https://mooseframework.inl.gov
Public Types | Public Member Functions | Protected Member Functions | Protected Attributes | Private Member Functions | Private Attributes | List of all members
MorrisReporterContext< DataType > Class Template Reference

#include <MorrisReporter.h>

Inheritance diagram for MorrisReporterContext< DataType >:
[legend]

Public Types

enum  AutoOperation
 

Public Member Functions

 MorrisReporterContext (const libMesh::ParallelObject &other, const MooseObject &producer, ReporterState< MorrisState< DataType >> &state, Sampler &sampler, const std::vector< DataType > &data)
 
 MorrisReporterContext (const libMesh::ParallelObject &other, const MooseObject &producer, ReporterState< MorrisState< DataType >> &state, Sampler &sampler, const std::vector< DataType > &data, const MooseEnum &ci_method, const std::vector< Real > &ci_levels, unsigned int ci_replicates, unsigned int ci_seed)
 
virtual void finalize () override
 
virtual std::string type () const override
 
virtual void declareClone (ReporterData &r_data, const ReporterName &r_name, const ReporterMode &mode, const MooseObject &producer) const override
 
virtual void declareVectorClone (ReporterData &r_data, const ReporterName &r_name, const ReporterMode &mode, const MooseObject &producer) const override
 
virtual void resize (dof_id_type local_size) final
 
virtual void clear () final
 
virtual void vectorSum () final
 
virtual std::string contextType () const override
 
const ReporterNamename () const override final
 
const ReporterState< MorrisState< DataType > > & state () const
 
virtual void transfer (ReporterData &r_data, const ReporterName &r_name, unsigned int time_index=0) const override
 
virtual void transferToVector (ReporterData &r_data, const ReporterName &r_name, dof_id_type index, unsigned int time_index=0) const override
 
virtual void transferFromVector (ReporterData &r_data, const ReporterName &r_name, dof_id_type index, unsigned int time_index=0) const override
 
void init (const ReporterMode &mode)
 
const MooseObjectgetProducer () const
 
const ReporterProducerEnumgetProducerModeEnum () const
 
const Parallel::Communicator & comm () const
 
processor_id_type n_processors () const
 
processor_id_type processor_id () const
 

Protected Member Functions

void broadcast ()
 
virtual void storeInfo (nlohmann::json &json) const override
 
virtual void store (nlohmann::json &json) const override
 
virtual void copyValuesBack () override
 
virtual bool restoreState () override
 
void requiresConsumerModes (const ReporterStateBase &state, const std::set< ReporterMode > &modes) const
 

Protected Attributes

ReporterState< MorrisState< DataType > > & _state
 
const MooseObject_producer
 
ReporterProducerEnum _producer_enum
 
const Parallel::Communicator & _communicator
 

Private Member Functions

std::vector< DataType > computeElementaryEffects (const RealEigenMatrix &x, const std::vector< DataType > &y) const
 Function for computing elementary effects for a single set of trajectories This is meant to be specialized for different data types. More...
 
template<>
std::vector< RealcomputeElementaryEffects (const RealEigenMatrix &x, const std::vector< Real > &y) const
 
template<>
std::vector< std::vector< Real > > computeElementaryEffects (const RealEigenMatrix &x, const std::vector< std::vector< Real >> &y) const
 

Private Attributes

Sampler_sampler
 Morris sampler (don't need any specific functions, but should be this type) More...
 
const std::vector< DataType > & _data
 Data used for the statistic calculation. More...
 
std::unique_ptr< StochasticTools::Calculator< std::vector< DataType >, DataType > > _mu_calc
 Storage for the Calculator object for the desired stat, this is created in constructor. More...
 
std::unique_ptr< StochasticTools::Calculator< std::vector< DataType >, DataType > > _mustar_calc
 
std::unique_ptr< StochasticTools::Calculator< std::vector< DataType >, DataType > > _sig_calc
 
std::unique_ptr< StochasticTools::BootstrapCalculator< std::vector< DataType >, DataType > > _ci_mu_calc
 Storage for the BootstrapCalculator for the desired confidence interval calculations (optional) More...
 
std::unique_ptr< StochasticTools::BootstrapCalculator< std::vector< DataType >, DataType > > _ci_mustar_calc
 
std::unique_ptr< StochasticTools::BootstrapCalculator< std::vector< DataType >, DataType > > _ci_sig_calc
 

Detailed Description

template<typename DataType>
class MorrisReporterContext< DataType >

Definition at line 57 of file MorrisReporter.h.

Constructor & Destructor Documentation

◆ MorrisReporterContext() [1/2]

template<typename DataType >
MorrisReporterContext< DataType >::MorrisReporterContext ( const libMesh::ParallelObject other,
const MooseObject producer,
ReporterState< MorrisState< DataType >> &  state,
Sampler sampler,
const std::vector< DataType > &  data 
)

Definition at line 147 of file MorrisReporter.C.

153  _sampler(sampler),
154  _data(data)
155 {
156  MultiMooseEnum items("mean meanabs stddev", "mean meanabs stddev", true);
157  _mu_calc = StochasticTools::makeCalculator<std::vector<DataType>, DataType>(items[0], other);
158  _mustar_calc = StochasticTools::makeCalculator<std::vector<DataType>, DataType>(items[1], other);
159  _sig_calc = StochasticTools::makeCalculator<std::vector<DataType>, DataType>(items[2], other);
160 
161  // Initialize state
162  auto & mu = this->_state.value()["mu"].first;
163  mu.resize(_sampler.getNumberOfCols());
164  auto & mu_star = this->_state.value()["mu_star"].first;
165  mu_star.resize(_sampler.getNumberOfCols());
166  auto & sig = this->_state.value()["sigma"].first;
167  sig.resize(_sampler.getNumberOfCols());
168 }
const ReporterState< MorrisState< DataType > > & state() const
ReporterState< MorrisState< DataType > > & _state
Sampler & _sampler
Morris sampler (don&#39;t need any specific functions, but should be this type)
static const std::string mu
Definition: NS.h:123
std::unique_ptr< StochasticTools::Calculator< std::vector< DataType >, DataType > > _mu_calc
Storage for the Calculator object for the desired stat, this is created in constructor.
const std::vector< DataType > & _data
Data used for the statistic calculation.
std::unique_ptr< StochasticTools::Calculator< std::vector< DataType >, DataType > > _sig_calc
T & value(const std::size_t time_index=0)
dof_id_type getNumberOfCols() const
std::unique_ptr< StochasticTools::Calculator< std::vector< DataType >, DataType > > _mustar_calc

◆ MorrisReporterContext() [2/2]

template<typename DataType >
MorrisReporterContext< DataType >::MorrisReporterContext ( const libMesh::ParallelObject other,
const MooseObject producer,
ReporterState< MorrisState< DataType >> &  state,
Sampler sampler,
const std::vector< DataType > &  data,
const MooseEnum ci_method,
const std::vector< Real > &  ci_levels,
unsigned int  ci_replicates,
unsigned int  ci_seed 
)

Definition at line 171 of file MorrisReporter.C.

180  : MorrisReporterContext<DataType>(other, producer, state, sampler, data)
181 {
182  _ci_mu_calc = StochasticTools::makeBootstrapCalculator<std::vector<DataType>, DataType>(
183  ci_method, other, ci_levels, ci_replicates, ci_seed, *_mu_calc);
184  _ci_mustar_calc = StochasticTools::makeBootstrapCalculator<std::vector<DataType>, DataType>(
185  ci_method, other, ci_levels, ci_replicates, ci_seed, *_mustar_calc);
186  _ci_sig_calc = StochasticTools::makeBootstrapCalculator<std::vector<DataType>, DataType>(
187  ci_method, other, ci_levels, ci_replicates, ci_seed, *_sig_calc);
188 
189  // Initialize state
190  auto & mu = this->_state.value()["mu"].second;
191  mu.resize(_sampler.getNumberOfCols());
192  auto & mu_star = this->_state.value()["mu_star"].second;
193  mu_star.resize(_sampler.getNumberOfCols());
194  auto & sig = this->_state.value()["sigma"].second;
195  sig.resize(_sampler.getNumberOfCols());
196 }
const ReporterState< MorrisState< DataType > > & state() const
ReporterState< MorrisState< DataType > > & _state
std::unique_ptr< StochasticTools::BootstrapCalculator< std::vector< DataType >, DataType > > _ci_mustar_calc
Sampler & _sampler
Morris sampler (don&#39;t need any specific functions, but should be this type)
static const std::string mu
Definition: NS.h:123
std::unique_ptr< StochasticTools::Calculator< std::vector< DataType >, DataType > > _mu_calc
Storage for the Calculator object for the desired stat, this is created in constructor.
std::unique_ptr< StochasticTools::Calculator< std::vector< DataType >, DataType > > _sig_calc
T & value(const std::size_t time_index=0)
dof_id_type getNumberOfCols() const
std::unique_ptr< StochasticTools::BootstrapCalculator< std::vector< DataType >, DataType > > _ci_sig_calc
std::unique_ptr< StochasticTools::BootstrapCalculator< std::vector< DataType >, DataType > > _ci_mu_calc
Storage for the BootstrapCalculator for the desired confidence interval calculations (optional) ...
std::unique_ptr< StochasticTools::Calculator< std::vector< DataType >, DataType > > _mustar_calc

Member Function Documentation

◆ computeElementaryEffects() [1/3]

template<typename DataType>
std::vector<DataType> MorrisReporterContext< DataType >::computeElementaryEffects ( const RealEigenMatrix x,
const std::vector< DataType > &  y 
) const
private

Function for computing elementary effects for a single set of trajectories This is meant to be specialized for different data types.

◆ computeElementaryEffects() [2/3]

template<>
std::vector< Real > MorrisReporterContext< Real >::computeElementaryEffects ( const RealEigenMatrix x,
const std::vector< Real > &  y 
) const
private

Definition at line 263 of file MorrisReporter.C.

265 {
266  const auto k = y.size() - 1;
267  const RealEigenMatrix dx = x.bottomRows(k) - x.topRows(k);
268  RealEigenVector dy(k);
269  for (unsigned int j = 0; j < k; ++j)
270  dy(j) = y[j + 1] - y[j];
271 
272  const RealEigenVector u = dx.fullPivLu().solve(dy);
273  return std::vector<Real>(u.data(), u.data() + u.size());
274 }
const std::vector< double > y
const std::vector< double > x
Eigen::Matrix< Real, Eigen::Dynamic, Eigen::Dynamic > RealEigenMatrix
static const std::complex< double > j(0, 1)
Complex number "j" (also known as "i")
Eigen::Matrix< Real, Eigen::Dynamic, 1 > RealEigenVector
static const std::string k
Definition: NS.h:130

◆ computeElementaryEffects() [3/3]

template<>
std::vector< std::vector< Real > > MorrisReporterContext< std::vector< Real > >::computeElementaryEffects ( const RealEigenMatrix x,
const std::vector< std::vector< Real >> &  y 
) const
private

Definition at line 278 of file MorrisReporter.C.

280 {
281  const auto k = y.size() - 1;
282  const RealEigenMatrix dx = x.bottomRows(k) - x.topRows(k);
283  const auto solver = dx.fullPivLu();
284  RealEigenVector dy(k);
285  std::vector<std::vector<Real>> ee(k, std::vector<Real>(y[0].size()));
286  for (unsigned int i = 0; i < y[0].size(); ++i)
287  {
288  for (unsigned int j = 0; j < k; ++j)
289  dy(j) = y[j + 1][i] - y[j][i];
290  const RealEigenVector u = solver.solve(dy);
291  for (unsigned int j = 0; j < k; ++j)
292  ee[j][i] = u(j);
293  }
294  return ee;
295 }
const std::vector< double > y
const std::vector< double > x
Eigen::Matrix< Real, Eigen::Dynamic, Eigen::Dynamic > RealEigenMatrix
static const std::complex< double > j(0, 1)
Complex number "j" (also known as "i")
Eigen::Matrix< Real, Eigen::Dynamic, 1 > RealEigenVector
static const std::string k
Definition: NS.h:130

◆ finalize()

template<typename DataType >
void MorrisReporterContext< DataType >::finalize ( )
overridevirtual

Reimplemented from ReporterGeneralContext< MorrisState< DataType > >.

Definition at line 200 of file MorrisReporter.C.

201 {
202  dof_id_type offset;
203  if (_data.size() == _sampler.getNumberOfRows())
204  offset = _sampler.getLocalRowBegin();
205  else if (_data.size() == _sampler.getNumberOfLocalRows())
206  offset = 0;
207  else
208  mooseError("Data size inconsistency. Expected data vector to have ",
210  " (local) or ",
212  " (global) elements, but actually has ",
213  _data.size(),
214  " elements. Are you using the same sampler?");
215 
216  const dof_id_type nc = _sampler.getNumberOfCols(); // convenience
217  if (_sampler.getNumberOfLocalRows() % (nc + 1) > 0)
218  mooseError(
219  "Sampler does not have an appropriate number of rows. Are you using a Morris sampler?");
220 
221  std::vector<std::vector<DataType>> elem_effects(
222  nc, std::vector<DataType>(_sampler.getNumberOfLocalRows() / (nc + 1)));
223  RealEigenMatrix x(nc + 1, nc);
224  std::vector<DataType> y(nc + 1);
225 
226  for (dof_id_type p = 0; p < _sampler.getNumberOfLocalRows(); ++p)
227  {
228  dof_id_type traj = p / (nc + 1);
229  dof_id_type tind = p % (nc + 1);
230  const std::vector<Real> row = _sampler.getNextLocalRow();
231  for (unsigned int k = 0; k < nc; ++k)
232  x(tind, k) = row[k];
233  y[tind] = _data[p + offset];
234 
235  if (tind == nc)
236  {
237  const std::vector<DataType> ee = computeElementaryEffects(x, y);
238  for (unsigned int k = 0; k < nc; ++k)
239  elem_effects[k][traj] = ee[k];
240  }
241  }
242 
243  auto & mu = this->_state.value()["mu"];
244  auto & mustar = this->_state.value()["mu_star"];
245  auto & sig = this->_state.value()["sigma"];
246  for (unsigned int k = 0; k < nc; ++k)
247  {
248  mu.first[k] = _mu_calc->compute(elem_effects[k], true);
249  mustar.first[k] = _mustar_calc->compute(elem_effects[k], true);
250  sig.first[k] = _sig_calc->compute(elem_effects[k], true);
251 
252  if (_ci_mu_calc)
253  mu.second[k] = _ci_mu_calc->compute(elem_effects[k], true);
254  if (_ci_mustar_calc)
255  mustar.second[k] = _ci_mustar_calc->compute(elem_effects[k], true);
256  if (_ci_sig_calc)
257  sig.second[k] = _ci_sig_calc->compute(elem_effects[k], true);
258  }
259 }
void mooseError(Args &&... args)
ReporterState< MorrisState< DataType > > & _state
std::vector< Real > getNextLocalRow()
dof_id_type getLocalRowBegin() const
std::unique_ptr< StochasticTools::BootstrapCalculator< std::vector< DataType >, DataType > > _ci_mustar_calc
const std::vector< double > y
dof_id_type getNumberOfLocalRows() const
std::vector< DataType > computeElementaryEffects(const RealEigenMatrix &x, const std::vector< DataType > &y) const
Function for computing elementary effects for a single set of trajectories This is meant to be specia...
const std::vector< double > x
Sampler & _sampler
Morris sampler (don&#39;t need any specific functions, but should be this type)
static const std::string mu
Definition: NS.h:123
Eigen::Matrix< Real, Eigen::Dynamic, Eigen::Dynamic > RealEigenMatrix
dof_id_type getNumberOfRows() const
std::unique_ptr< StochasticTools::Calculator< std::vector< DataType >, DataType > > _mu_calc
Storage for the Calculator object for the desired stat, this is created in constructor.
const std::vector< DataType > & _data
Data used for the statistic calculation.
std::unique_ptr< StochasticTools::Calculator< std::vector< DataType >, DataType > > _sig_calc
T & value(const std::size_t time_index=0)
static const std::string k
Definition: NS.h:130
dof_id_type getNumberOfCols() const
uint8_t dof_id_type
std::unique_ptr< StochasticTools::BootstrapCalculator< std::vector< DataType >, DataType > > _ci_sig_calc
std::unique_ptr< StochasticTools::BootstrapCalculator< std::vector< DataType >, DataType > > _ci_mu_calc
Storage for the BootstrapCalculator for the desired confidence interval calculations (optional) ...
std::unique_ptr< StochasticTools::Calculator< std::vector< DataType >, DataType > > _mustar_calc

◆ type()

template<typename DataType>
virtual std::string MorrisReporterContext< DataType >::type ( ) const
inlineoverridevirtual

Reimplemented from ReporterGeneralContext< MorrisState< DataType > >.

Definition at line 77 of file MorrisReporter.h.

78  {
79  return "MorrisSensitivity<" + MooseUtils::prettyCppType<DataType>() + ">";
80  }

Member Data Documentation

◆ _ci_mu_calc

template<typename DataType>
std::unique_ptr<StochasticTools::BootstrapCalculator<std::vector<DataType>, DataType> > MorrisReporterContext< DataType >::_ci_mu_calc
private

Storage for the BootstrapCalculator for the desired confidence interval calculations (optional)

Definition at line 103 of file MorrisReporter.h.

Referenced by MorrisReporterContext< DataType >::MorrisReporterContext().

◆ _ci_mustar_calc

template<typename DataType>
std::unique_ptr<StochasticTools::BootstrapCalculator<std::vector<DataType>, DataType> > MorrisReporterContext< DataType >::_ci_mustar_calc
private

◆ _ci_sig_calc

template<typename DataType>
std::unique_ptr<StochasticTools::BootstrapCalculator<std::vector<DataType>, DataType> > MorrisReporterContext< DataType >::_ci_sig_calc
private

◆ _data

template<typename DataType>
const std::vector<DataType>& MorrisReporterContext< DataType >::_data
private

Data used for the statistic calculation.

Definition at line 94 of file MorrisReporter.h.

◆ _mu_calc

template<typename DataType>
std::unique_ptr<StochasticTools::Calculator<std::vector<DataType>, DataType> > MorrisReporterContext< DataType >::_mu_calc
private

Storage for the Calculator object for the desired stat, this is created in constructor.

Definition at line 97 of file MorrisReporter.h.

Referenced by MorrisReporterContext< DataType >::MorrisReporterContext().

◆ _mustar_calc

template<typename DataType>
std::unique_ptr<StochasticTools::Calculator<std::vector<DataType>, DataType> > MorrisReporterContext< DataType >::_mustar_calc
private

◆ _sampler

template<typename DataType>
Sampler& MorrisReporterContext< DataType >::_sampler
private

Morris sampler (don't need any specific functions, but should be this type)

Definition at line 91 of file MorrisReporter.h.

Referenced by MorrisReporterContext< DataType >::MorrisReporterContext().

◆ _sig_calc

template<typename DataType>
std::unique_ptr<StochasticTools::Calculator<std::vector<DataType>, DataType> > MorrisReporterContext< DataType >::_sig_calc
private

The documentation for this class was generated from the following files: