Stokhos Package Browser (Single Doxygen Collection)  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Stokhos_SGModelEvaluator.cpp
Go to the documentation of this file.
1 // @HEADER
2 // *****************************************************************************
3 // Stokhos Package
4 //
5 // Copyright 2009 NTESS and the Stokhos contributors.
6 // SPDX-License-Identifier: BSD-3-Clause
7 // *****************************************************************************
8 // @HEADER
9 
11 
12 #include <algorithm>
13 #include "Teuchos_Assert.hpp"
14 #include "EpetraExt_BlockUtility.h"
15 #include "EpetraExt_BlockMultiVector.h"
20 #include "Epetra_LocalMap.h"
21 #include "Epetra_Export.h"
22 #include "Epetra_Import.h"
23 
26  const Teuchos::RCP<const Stokhos::OrthogPolyBasis<int,double> >& sg_basis_,
27  const Teuchos::RCP<const Stokhos::Quadrature<int,double> >& sg_quad_,
29  const Teuchos::RCP<const Stokhos::ParallelData>& sg_parallel_data_,
31  bool scaleOP_)
32  : me(me_),
33  sg_basis(sg_basis_),
34  sg_quad(sg_quad_),
35  sg_exp(sg_exp_),
36  params(params_),
37  num_sg_blocks(sg_basis->size()),
38  num_W_blocks(sg_basis->size()),
39  num_p_blocks(sg_basis->size()),
40  supports_x(false),
41  x_map(me->get_x_map()),
42  f_map(me->get_f_map()),
43  sg_parallel_data(sg_parallel_data_),
44  sg_comm(sg_parallel_data->getMultiComm()),
45  epetraCijk(sg_parallel_data->getEpetraCijk()),
46  serialCijk(),
47  sg_x_map(),
48  sg_overlapped_x_map(),
49  sg_f_map(),
50  sg_overlapped_f_map(),
51  sg_overlapped_x_importer(),
52  sg_overlapped_f_exporter(),
53  num_p(0),
54  num_p_sg(0),
55  sg_p_index_map(),
56  sg_p_map(),
57  sg_p_names(),
58  num_g(0),
59  num_g_sg(0),
60  sg_g_index_map(),
61  sg_g_map(),
62  x_dot_sg_blocks(),
63  x_sg_blocks(),
64  f_sg_blocks(),
65  W_sg_blocks(),
66  sg_x_init(),
67  sg_p_init(),
68  eval_W_with_f(false),
69  scaleOP(scaleOP_)
70 {
71  if (x_map != Teuchos::null)
72  supports_x = true;
73 
76  static_cast<int>(num_sg_blocks), 0, *(sg_parallel_data->getStochasticComm())));
79  stoch_row_map = epetraCijk->getStochasticRowMap();
80 
81  if (supports_x) {
82 
83  // Create block SG x and f maps
84  sg_x_map =
85  Teuchos::rcp(EpetraExt::BlockUtility::GenerateBlockMap(
86  *x_map, *stoch_row_map, *sg_comm));
88  Teuchos::rcp(EpetraExt::BlockUtility::GenerateBlockMap(
90 
91  sg_f_map =
92  Teuchos::rcp(EpetraExt::BlockUtility::GenerateBlockMap(
93  *f_map, *stoch_row_map, *sg_comm));
95  Teuchos::rcp(EpetraExt::BlockUtility::GenerateBlockMap(
97 
98  // Create importer/exporter from/to overlapped distribution
103 
104  // Create vector blocks
108 
109  // Create default sg_x_init
111  if (sg_x_init->myGID(0))
112  (*sg_x_init)[0] = *(me->get_x_init());
113 
114  // Preconditioner needs an x
116 
117  // Determine W expansion type
118  std::string W_expansion_type =
119  params->get("Jacobian Expansion Type", "Full");
120  if (W_expansion_type == "Linear")
122  else
124  Teuchos::RCP<Epetra_BlockMap> W_overlap_map =
126  static_cast<int>(num_W_blocks), 0,
127  *(sg_parallel_data->getStochasticComm())));
128  W_sg_blocks =
130  sg_basis, W_overlap_map, x_map, f_map, sg_f_map, sg_comm));
131  for (unsigned int i=0; i<num_W_blocks; i++)
132  W_sg_blocks->setCoeffPtr(i, me->create_W());
133 
134  eval_W_with_f = params->get("Evaluate W with F", false);
135 
137  Teuchos::rcp(&(params->sublist("SG Preconditioner")), false);
140  }
141 
142  // Parameters -- The idea here is to add new parameter vectors
143  // for the stochastic Galerkin components of the parameters
144 
145  InArgs me_inargs = me->createInArgs();
146  OutArgs me_outargs = me->createOutArgs();
147  num_p = me_inargs.Np();
148 
149  // Get the p_sg's supported and build index map
150  for (int i=0; i<num_p; i++) {
151  if (me_inargs.supports(IN_ARG_p_sg, i))
153  }
155 
156  sg_p_map.resize(num_p_sg);
157  sg_p_names.resize(num_p_sg);
158  sg_p_init.resize(num_p_sg);
159 
160  // Determine parameter expansion type
161  std::string p_expansion_type =
162  params->get("Parameter Expansion Type", "Full");
163  if (p_expansion_type == "Linear")
165  else
167 
168  // Create parameter maps, names, and initial values
171  static_cast<int>(num_p_blocks), 0,
172  *(sg_parallel_data->getStochasticComm())));
173  for (int i=0; i<num_p_sg; i++) {
174  Teuchos::RCP<const Epetra_Map> p_map = me->get_p_map(sg_p_index_map[i]);
175  sg_p_map[i] =
176  Teuchos::rcp(EpetraExt::BlockUtility::GenerateBlockMap(
177  *p_map, *overlapped_stoch_p_map, *sg_comm));
178 
180  me->get_p_names(sg_p_index_map[i]);
181  if (p_names != Teuchos::null) {
182  sg_p_names[i] =
184  for (int j=0; j<p_names->size(); j++) {
185  std::stringstream ss;
186  ss << (*p_names)[j] << " -- SG Coefficient " << i;
187  (*sg_p_names[i])[j] = ss.str();
188  }
189  }
190 
191  // Create default sg_p_init with an expansion equal to the mean parameter
193  sg_p_init[i]->init(0.0);
194  (*sg_p_init[i])[0] = *(me->get_p_init(sg_p_index_map[i]));
195  }
196 
197  // Responses -- The idea here is to add new parameter vectors
198  // for the stochastic Galerkin components of the respones
199  num_g = me_outargs.Ng();
200 
201  // Get the g_sg's supported and build index map
202  for (int i=0; i<num_g; i++) {
203  if (me_outargs.supports(OUT_ARG_g_sg, i))
205  }
207 
208  sg_g_map.resize(num_g_sg);
209 
210  // Create response maps
211  for (int i=0; i<num_g_sg; i++) {
212  Teuchos::RCP<const Epetra_Map> g_map = me->get_g_map(sg_g_index_map[i]);
213  sg_g_map[i] =
214  Teuchos::rcp(EpetraExt::BlockUtility::GenerateBlockMap(
215  *g_map, *overlapped_stoch_row_map, *sg_comm));
216  }
217 
218  // We don't support parallel for dgdx yet, so build a new EpetraCijk
219  if (supports_x) {
220  serialCijk =
222  epetraCijk->getCijk(),
223  sg_comm,
225  }
226 
227 }
228 
229 // Overridden from EpetraExt::ModelEvaluator
230 
233 {
234  return sg_x_map;
235 }
236 
239 {
240  return sg_f_map;
241 }
242 
245 {
246  TEUCHOS_TEST_FOR_EXCEPTION(l < 0 || l >= num_p + num_p_sg, std::logic_error,
247  "Error! Invalid p map index " << l);
248  if (l < num_p)
249  return me->get_p_map(l);
250  else
251  return sg_p_map[l-num_p];
252 
253  return Teuchos::null;
254 }
255 
258 {
259  TEUCHOS_TEST_FOR_EXCEPTION(j < 0 || j >= num_g_sg, std::logic_error,
260  "Error! Invalid g map index " << j);
261  return sg_g_map[j];
262 }
263 
266 {
267  TEUCHOS_TEST_FOR_EXCEPTION(l < 0 || l >= num_p + num_p_sg, std::logic_error,
268  "Error! Invalid p map index " << l);
269  if (l < num_p)
270  return me->get_p_names(l);
271  else
272  return sg_p_names[l-num_p];
273 
274  return Teuchos::null;
275 }
276 
279 {
280  return sg_x_init->getBlockVector();
281 }
282 
285 {
286  TEUCHOS_TEST_FOR_EXCEPTION(l < 0 || l >= num_p + num_p_sg, std::logic_error,
287  "Error! Invalid p map index " << l);
288  if (l < num_p)
289  return me->get_p_init(l);
290  else
291  return sg_p_init[l-num_p]->getBlockVector();
292 
293  return Teuchos::null;
294 }
295 
298 {
299  if (supports_x) {
301  Teuchos::rcp(&(params->sublist("SG Operator")), false);
303  if (sgOpParams->get("Operator Method", "Matrix Free") ==
304  "Fully Assembled") {
305  W_crs = Teuchos::rcp_dynamic_cast<Epetra_CrsMatrix>(me->create_W(), true);
307  Teuchos::rcp(&(W_crs->Graph()), false);
308  sgOpParams->set< Teuchos::RCP<const Epetra_CrsGraph> >("Base Graph",
309  W_graph);
310  }
311  Stokhos::SGOperatorFactory sg_op_factory(sgOpParams);
312  my_W = sg_op_factory.build(sg_comm, sg_basis, epetraCijk, x_map, f_map,
313  sg_x_map, sg_f_map);
314  my_W->setupOperator(W_sg_blocks);
315 
316  return my_W;
317  }
318 
319  return Teuchos::null;
320 }
321 
324 {
325  if (supports_x) {
326 
328  sg_prec_factory->build(sg_comm, sg_basis, epetraCijk, x_map, sg_x_map);
329  return Teuchos::rcp(new EpetraExt::ModelEvaluator::Preconditioner(precOp,
330  true));
331  }
332 
333  return Teuchos::null;
334 }
335 
338 {
339  TEUCHOS_TEST_FOR_EXCEPTION(j < 0 || j >= num_g_sg || !supports_x,
340  std::logic_error,
341  "Error: dg/dx index " << j << " is not supported!");
342 
343  int jj = sg_g_index_map[j];
344  Teuchos::RCP<const Epetra_Map> g_map = me->get_g_map(jj);
346  OutArgs me_outargs = me->createOutArgs();
347  DerivativeSupport ds = me_outargs.supports(OUT_ARG_DgDx_sg, jj);
348  if (ds.supports(DERIV_LINEAR_OP)) {
349  sg_blocks =
351  sg_basis, overlapped_stoch_row_map, x_map,
352  g_map, sg_g_map[j], sg_comm));
353  for (unsigned int i=0; i<num_sg_blocks; i++)
354  sg_blocks->setCoeffPtr(i, me->create_DgDx_op(jj));
355  }
356  else if (ds.supports(DERIV_MV_BY_COL)) {
359  sg_basis, overlapped_stoch_row_map, g_map, sg_g_map[j],
360  sg_comm, x_map->NumMyElements()));
361  sg_blocks =
363  sg_mv_blocks, false));
364  }
365  else if (ds.supports(DERIV_TRANS_MV_BY_ROW)) {
368  sg_basis, overlapped_stoch_row_map, x_map, sg_x_map,
369  sg_comm, g_map->NumMyElements()));
370  sg_blocks =
372  sg_mv_blocks, true));
373  }
374  else
375  TEUCHOS_TEST_FOR_EXCEPTION(true, std::logic_error,
376  "Error! me_outargs.supports(OUT_ARG_DgDx_sg, " << jj
377  << ").none() is true!");
378 
383  sg_comm, sg_basis, serialCijk, x_map, g_map,
384  sg_x_map, sg_g_map[j], pl));
385  dgdx_sg->setupOperator(sg_blocks);
386  return dgdx_sg;
387 }
388 
391 {
392  TEUCHOS_TEST_FOR_EXCEPTION(j < 0 || j >= num_g_sg || !supports_x,
393  std::logic_error,
394  "Error: dg/dx_dot index " << j << " is not supported!");
395 
396  int jj = sg_g_index_map[j];
397  Teuchos::RCP<const Epetra_Map> g_map = me->get_g_map(jj);
399  OutArgs me_outargs = me->createOutArgs();
400  DerivativeSupport ds = me_outargs.supports(OUT_ARG_DgDx_dot_sg, jj);
401  if (ds.supports(DERIV_LINEAR_OP)) {
402  sg_blocks =
404  sg_basis, overlapped_stoch_row_map, x_map,
405  g_map, sg_g_map[j], sg_comm));
406  for (unsigned int i=0; i<num_sg_blocks; i++)
407  sg_blocks->setCoeffPtr(i, me->create_DgDx_dot_op(jj));
408  }
409  else if (ds.supports(DERIV_MV_BY_COL)) {
412  sg_basis, overlapped_stoch_row_map, g_map, sg_g_map[j],
413  sg_comm, x_map->NumMyElements()));
414  sg_blocks =
416  sg_mv_blocks, false));
417  }
418  else if (ds.supports(DERIV_TRANS_MV_BY_ROW)) {
421  sg_basis, overlapped_stoch_row_map, x_map, sg_x_map,
422  sg_comm, g_map->NumMyElements()));
423  sg_blocks =
425  sg_mv_blocks, true));
426  }
427  else
428  TEUCHOS_TEST_FOR_EXCEPTION(true, std::logic_error,
429  "Error! me_outargs.supports(OUT_ARG_DgDx_dot_sg, "
430  << jj << ").none() is true!");
431 
436  sg_comm, sg_basis, serialCijk, x_map, g_map,
437  sg_x_map, sg_g_map[j], pl));
438  dgdx_dot_sg->setupOperator(sg_blocks);
439  return dgdx_dot_sg;
440 }
441 
444 {
446  j < 0 || j >= num_g_sg || i < 0 || i >= num_p+num_p_sg,
447  std::logic_error,
448  "Error: dg/dp index " << j << "," << i << " is not supported!");
449 
450  int jj = sg_g_index_map[j];
451  Teuchos::RCP<const Epetra_Map> g_map = me->get_g_map(jj);
452  OutArgs me_outargs = me->createOutArgs();
453  if (i < num_p) {
454  if (me_outargs.supports(OUT_ARG_DgDp_sg,jj,i).supports(DERIV_LINEAR_OP)) {
457  sg_basis, overlapped_stoch_row_map,
458  me->get_p_map(i), me->get_g_map(jj),
459  sg_g_map[j], sg_comm));
460  for (unsigned int l=0; l<num_sg_blocks; l++)
461  sg_blocks->setCoeffPtr(l, me->create_DgDp_op(jj,i));
462  return sg_blocks;
463  }
464  else
466  true, std::logic_error,
467  "Error: Underlying model evaluator must support DERIV_LINER_OP " <<
468  "to create operator for dg/dp index " << j << "," << i << "!");
469  }
470  else {
471  int ii = sg_p_index_map[i-num_p];
472  Teuchos::RCP<const Epetra_Map> p_map = me->get_p_map(ii);
474  DerivativeSupport ds = me_outargs.supports(OUT_ARG_DgDp_sg,jj,ii);
475  if (ds.supports(DERIV_LINEAR_OP)) {
476  sg_blocks =
478  sg_basis, overlapped_stoch_row_map,
479  p_map, g_map, sg_g_map[j], sg_comm));
480  for (unsigned int l=0; l<num_sg_blocks; l++)
481  sg_blocks->setCoeffPtr(l, me->create_DgDp_op(jj,ii));
482  }
483  else if (ds.supports(DERIV_MV_BY_COL)) {
486  sg_basis, overlapped_stoch_row_map, g_map, sg_g_map[j],
487  sg_comm, p_map->NumMyElements()));
488  sg_blocks =
490  sg_mv_blocks, false));
491  }
492  else if (ds.supports(DERIV_TRANS_MV_BY_ROW)) {
495  sg_basis, overlapped_stoch_row_map, p_map,
496  sg_p_map[i-num_p],
497  sg_comm, g_map->NumMyElements()));
498  sg_blocks =
500  sg_mv_blocks, true));
501  }
502  else
503  TEUCHOS_TEST_FOR_EXCEPTION(true, std::logic_error,
504  "Error! me_outargs.supports(OUT_ARG_DgDp_sg, " << jj
505  << "," << ii << ").none() is true!");
506 
511  sg_comm, sg_basis, serialCijk,
512  p_map, g_map, sg_p_map[i-num_p], sg_g_map[j], pl));
513  dgdp_sg->setupOperator(sg_blocks);
514  return dgdp_sg;
515  }
516 
517  return Teuchos::null;
518 }
519 
522 {
523  TEUCHOS_TEST_FOR_EXCEPTION(i < 0 || i >= num_p+num_p_sg,
524  std::logic_error,
525  "Error: df/dp index " << i << " is not supported!");
526 
527  OutArgs me_outargs = me->createOutArgs();
528  if (i < num_p) {
529  if (me_outargs.supports(OUT_ARG_DfDp_sg,i).supports(DERIV_LINEAR_OP)) {
532  sg_basis, overlapped_stoch_row_map,
533  me->get_p_map(i), f_map,
534  sg_overlapped_f_map, sg_comm));
535  for (unsigned int l=0; l<num_sg_blocks; l++)
536  sg_blocks->setCoeffPtr(l, me->create_DfDp_op(i));
537  return sg_blocks;
538  }
539  else
541  true, std::logic_error,
542  "Error: Underlying model evaluator must support DERIV_LINER_OP " <<
543  "to create operator for df/dp index " << i << "!");
544  }
545  else {
546  int ii = sg_p_index_map[i-num_p];
547  Teuchos::RCP<const Epetra_Map> p_map = me->get_p_map(ii);
549  DerivativeSupport ds = me_outargs.supports(OUT_ARG_DfDp_sg,ii);
550  if (ds.supports(DERIV_LINEAR_OP)) {
551  sg_blocks =
553  sg_basis, overlapped_stoch_row_map,
554  p_map, f_map, sg_overlapped_f_map, sg_comm));
555  for (unsigned int l=0; l<num_sg_blocks; l++)
556  sg_blocks->setCoeffPtr(l, me->create_DfDp_op(ii));
557  }
558  else if (ds.supports(DERIV_MV_BY_COL)) {
561  sg_basis, overlapped_stoch_row_map, f_map, sg_f_map,
562  sg_comm, p_map->NumMyElements()));
563  sg_blocks =
565  sg_mv_blocks, false));
566  }
567  else if (ds.supports(DERIV_TRANS_MV_BY_ROW)) {
570  sg_basis, overlapped_stoch_row_map, p_map,
571  sg_p_map[i-num_p],
572  sg_comm, f_map->NumMyElements()));
573  sg_blocks =
575  sg_mv_blocks, true));
576  }
577  else
578  TEUCHOS_TEST_FOR_EXCEPTION(true, std::logic_error,
579  "Error! me_outargs.supports(OUT_ARG_DfDp_sg, " << ii
580  << ").none() is true!");
581 
586  sg_comm, sg_basis, epetraCijk,
587  p_map, f_map, sg_p_map[i-num_p], sg_f_map, pl));
588  dfdp_sg->setupOperator(sg_blocks);
589  return dfdp_sg;
590  }
591 
592  return Teuchos::null;
593 }
594 
595 EpetraExt::ModelEvaluator::InArgs
597 {
598  InArgsSetup inArgs;
599  InArgs me_inargs = me->createInArgs();
600 
601  inArgs.setModelEvalDescription(this->description());
602  inArgs.set_Np(num_p+num_p_sg);
603  inArgs.setSupports(IN_ARG_x_dot, me_inargs.supports(IN_ARG_x_dot_sg));
604  inArgs.setSupports(IN_ARG_x, me_inargs.supports(IN_ARG_x_sg));
605  inArgs.setSupports(IN_ARG_t, me_inargs.supports(IN_ARG_t));
606  inArgs.setSupports(IN_ARG_alpha, me_inargs.supports(IN_ARG_alpha));
607  inArgs.setSupports(IN_ARG_beta, me_inargs.supports(IN_ARG_beta));
608  inArgs.setSupports(IN_ARG_sg_basis, me_inargs.supports(IN_ARG_sg_basis));
609  inArgs.setSupports(IN_ARG_sg_quadrature,
610  me_inargs.supports(IN_ARG_sg_quadrature));
611  inArgs.setSupports(IN_ARG_sg_expansion,
612  me_inargs.supports(IN_ARG_sg_expansion));
613 
614  return inArgs;
615 }
616 
617 EpetraExt::ModelEvaluator::OutArgs
619 {
620  OutArgsSetup outArgs;
621  OutArgs me_outargs = me->createOutArgs();
622 
623  outArgs.setModelEvalDescription(this->description());
624  outArgs.set_Np_Ng(num_p+num_p_sg, num_g_sg);
625  outArgs.setSupports(OUT_ARG_f, me_outargs.supports(OUT_ARG_f_sg));
626  outArgs.setSupports(OUT_ARG_W, me_outargs.supports(OUT_ARG_W_sg));
627  outArgs.setSupports(
628  OUT_ARG_WPrec,
629  me_outargs.supports(OUT_ARG_W_sg) && sg_prec_factory->isPrecSupported());
630  for (int j=0; j<num_p; j++)
631  outArgs.setSupports(OUT_ARG_DfDp, j,
632  me_outargs.supports(OUT_ARG_DfDp_sg, j));
633  for (int j=0; j<num_p_sg; j++)
634  if (!me_outargs.supports(OUT_ARG_DfDp_sg, sg_p_index_map[j]).none())
635  outArgs.setSupports(OUT_ARG_DfDp, j+num_p, DERIV_LINEAR_OP);
636  for (int i=0; i<num_g_sg; i++) {
637  int ii = sg_g_index_map[i];
638  if (!me_outargs.supports(OUT_ARG_DgDx_dot_sg, ii).none())
639  outArgs.setSupports(OUT_ARG_DgDx_dot, i, DERIV_LINEAR_OP);
640  if (!me_outargs.supports(OUT_ARG_DgDx_sg, ii).none())
641  outArgs.setSupports(OUT_ARG_DgDx, i, DERIV_LINEAR_OP);
642  for (int j=0; j<num_p; j++)
643  outArgs.setSupports(OUT_ARG_DgDp, i, j,
644  me_outargs.supports(OUT_ARG_DgDp_sg, ii, j));
645  for (int j=0; j<num_p_sg; j++)
646  if (!me_outargs.supports(OUT_ARG_DgDp_sg, ii, sg_p_index_map[j]).none())
647  outArgs.setSupports(OUT_ARG_DgDp, i, j+num_p, DERIV_LINEAR_OP);
648  }
649 
650  return outArgs;
651 }
652 
653 void
655  const OutArgs& outArgs) const
656 {
657  // Get the input arguments
659  if (inArgs.supports(IN_ARG_x)) {
660  x = inArgs.get_x();
661  if (x != Teuchos::null)
662  *my_x = *x;
663  }
665  if (inArgs.supports(IN_ARG_x_dot))
666  x_dot = inArgs.get_x_dot();
667 
668  // Get the output arguments
669  EpetraExt::ModelEvaluator::Evaluation<Epetra_Vector> f_out;
670  if (outArgs.supports(OUT_ARG_f))
671  f_out = outArgs.get_f();
673  if (outArgs.supports(OUT_ARG_W))
674  W_out = outArgs.get_W();
676  if (outArgs.supports(OUT_ARG_WPrec))
677  WPrec_out = outArgs.get_WPrec();
678 
679  // Check if we are using the "matrix-free" method for W and we are
680  // computing a preconditioner.
681  bool eval_prec = (W_out == Teuchos::null && WPrec_out != Teuchos::null);
682 
683  // Here we are assuming a full W fill occurred previously which we can use
684  // for the preconditioner. Given the expense of computing the SG W blocks
685  // this saves significant computational cost
686  if (eval_prec) {
688  Teuchos::rcp_dynamic_cast<Stokhos::SGPreconditioner>(WPrec_out, true);
689  W_prec->setupPreconditioner(my_W, *my_x);
690 
691  // We can now quit unless a fill of f, g, or dg/dp was also requested
692  bool done = (f_out == Teuchos::null);
693  for (int i=0; i<outArgs.Ng(); i++) {
694  done = done && (outArgs.get_g(i) == Teuchos::null);
695  for (int j=0; j<outArgs.Np(); j++)
696  if (!outArgs.supports(OUT_ARG_DgDp, i, j).none())
697  done = done && (outArgs.get_DgDp(i,j).isEmpty());
698  }
699  if (done)
700  return;
701  }
702 
703  // Create underlying inargs
704  InArgs me_inargs = me->createInArgs();
705  if (x != Teuchos::null) {
706  x_sg_blocks->getBlockVector()->Import(*x, *sg_overlapped_x_importer,
707  Insert);
708  me_inargs.set_x_sg(x_sg_blocks);
709  }
710  if (x_dot != Teuchos::null) {
711  x_dot_sg_blocks->getBlockVector()->Import(*x_dot, *sg_overlapped_x_importer,
712  Insert);
713  me_inargs.set_x_dot_sg(x_dot_sg_blocks);
714  }
715  if (me_inargs.supports(IN_ARG_alpha))
716  me_inargs.set_alpha(inArgs.get_alpha());
717  if (me_inargs.supports(IN_ARG_beta))
718  me_inargs.set_beta(inArgs.get_beta());
719  if (me_inargs.supports(IN_ARG_t))
720  me_inargs.set_t(inArgs.get_t());
721  if (me_inargs.supports(IN_ARG_sg_basis)) {
722  if (inArgs.get_sg_basis() != Teuchos::null)
723  me_inargs.set_sg_basis(inArgs.get_sg_basis());
724  else
725  me_inargs.set_sg_basis(sg_basis);
726  }
727  if (me_inargs.supports(IN_ARG_sg_quadrature)) {
728  if (inArgs.get_sg_quadrature() != Teuchos::null)
729  me_inargs.set_sg_quadrature(inArgs.get_sg_quadrature());
730  else
731  me_inargs.set_sg_quadrature(sg_quad);
732  }
733  if (me_inargs.supports(IN_ARG_sg_expansion)) {
734  if (inArgs.get_sg_expansion() != Teuchos::null)
735  me_inargs.set_sg_expansion(inArgs.get_sg_expansion());
736  else
737  me_inargs.set_sg_expansion(sg_exp);
738  }
739 
740  // Pass parameters
741  for (int i=0; i<num_p; i++)
742  me_inargs.set_p(i, inArgs.get_p(i));
743  for (int i=0; i<num_p_sg; i++) {
744  Teuchos::RCP<const Epetra_Vector> p = inArgs.get_p(i+num_p);
745 
746  // We always need to pass in the SG parameters, so just use
747  // the initial parameters if it is NULL
748  if (p == Teuchos::null)
749  p = sg_p_init[i]->getBlockVector();
750 
751  // Convert block p to SG polynomial
753  create_p_sg(sg_p_index_map[i], View, p.get());
754  me_inargs.set_p_sg(sg_p_index_map[i], p_sg);
755  }
756 
757  // Create underlying outargs
758  OutArgs me_outargs = me->createOutArgs();
759 
760  // f
761  if (f_out != Teuchos::null) {
762  me_outargs.set_f_sg(f_sg_blocks);
763  if (eval_W_with_f)
764  me_outargs.set_W_sg(W_sg_blocks);
765  }
766 
767  // W
768  if (W_out != Teuchos::null && !eval_W_with_f && !eval_prec)
769  me_outargs.set_W_sg(W_sg_blocks);
770 
771  // df/dp -- deterministic p
772  for (int i=0; i<num_p; i++) {
773  if (!outArgs.supports(OUT_ARG_DfDp, i).none()) {
774  Derivative dfdp = outArgs.get_DfDp(i);
775  if (dfdp.getMultiVector() != Teuchos::null) {
777  if (dfdp.getMultiVectorOrientation() == DERIV_MV_BY_COL)
778  dfdp_sg =
780  sg_basis, overlapped_stoch_row_map,
781  me->get_f_map(), sg_overlapped_f_map, sg_comm,
782  me->get_p_map(i)->NumMyElements()));
783  else if (dfdp.getMultiVectorOrientation() == DERIV_TRANS_MV_BY_ROW)
784  dfdp_sg =
786  sg_basis, overlapped_stoch_row_map,
787  me->get_p_map(i), sg_comm,
788  me->get_f_map()->NumMyElements()));
789  me_outargs.set_DfDp_sg(
790  i, SGDerivative(dfdp_sg, dfdp.getMultiVectorOrientation()));
791  }
792  else if (dfdp.getLinearOp() != Teuchos::null) {
794  Teuchos::rcp_dynamic_cast<Stokhos::EpetraOperatorOrthogPoly>(dfdp.getLinearOp(), true);
795  me_outargs.set_DfDp_sg(i, SGDerivative(dfdp_sg));
796  }
797  }
798  }
799 
800  // dfdp -- stochastic p. Here we only support DERIV_LINEAR_OP
801  for (int i=0; i<num_p_sg; i++) {
802  if (!outArgs.supports(OUT_ARG_DfDp, i+num_p).none()) {
803  Derivative dfdp = outArgs.get_DfDp(i+num_p);
804  if (dfdp.getLinearOp() != Teuchos::null) {
806  Teuchos::rcp_dynamic_cast<Stokhos::MatrixFreeOperator>(dfdp.getLinearOp(), true);
808  dfdp_op->getSGPolynomial();
809  int ii = sg_p_index_map[i];
810  if (me_outargs.supports(OUT_ARG_DfDp_sg,ii).supports(DERIV_LINEAR_OP))
811  me_outargs.set_DfDp_sg(ii, SGDerivative(dfdp_op_sg));
812  else {
814  Teuchos::rcp_dynamic_cast<Stokhos::EpetraMultiVectorOperatorOrthogPoly>(dfdp_op_sg, true);
816  sg_mv_op->multiVectorOrthogPoly();
817  if (me_outargs.supports(OUT_ARG_DfDp_sg,ii).supports(DERIV_MV_BY_COL))
818  me_outargs.set_DfDp_sg(
819  ii, SGDerivative(dfdp_sg, DERIV_MV_BY_COL));
820  else
821  me_outargs.set_DfDp_sg(
822  ii, SGDerivative(dfdp_sg, DERIV_TRANS_MV_BY_ROW));
823  }
824  }
826  dfdp.getLinearOp() == Teuchos::null && dfdp.isEmpty() == false,
827  std::logic_error,
828  "Error! Stokhos::SGModelEvaluator::evalModel: " <<
829  "Operator form of df/dp(" << i+num_p << ") is required!");
830  }
831  }
832 
833  // Responses (g, dg/dx, dg/dp, ...)
834  for (int i=0; i<num_g_sg; i++) {
835  int ii = sg_g_index_map[i];
836 
837  // g
838  Teuchos::RCP<Epetra_Vector> g = outArgs.get_g(i);
839  if (g != Teuchos::null) {
841  create_g_sg(sg_g_index_map[i], View, g.get());
842  me_outargs.set_g_sg(ii, g_sg);
843  }
844 
845  // dg/dxdot
846  if (outArgs.supports(OUT_ARG_DgDx_dot, i).supports(DERIV_LINEAR_OP)) {
847  Derivative dgdx_dot = outArgs.get_DgDx_dot(i);
848  if (dgdx_dot.getLinearOp() != Teuchos::null) {
850  Teuchos::rcp_dynamic_cast<Stokhos::SGOperator>(
851  dgdx_dot.getLinearOp(), true);
853  op->getSGPolynomial();
854  if (me_outargs.supports(OUT_ARG_DgDx, ii).supports(DERIV_LINEAR_OP))
855  me_outargs.set_DgDx_dot_sg(ii, sg_blocks);
856  else {
858  Teuchos::rcp_dynamic_cast<Stokhos::EpetraMultiVectorOperatorOrthogPoly>(sg_blocks, true);
860  sg_mv_op->multiVectorOrthogPoly();
861  if (me_outargs.supports(OUT_ARG_DgDx_dot_sg, ii).supports(DERIV_MV_BY_COL))
862  me_outargs.set_DgDx_dot_sg(ii, SGDerivative(dgdx_dot_sg,
863  DERIV_MV_BY_COL));
864  else
865  me_outargs.set_DgDx_dot_sg(ii, SGDerivative(dgdx_dot_sg,
866  DERIV_TRANS_MV_BY_ROW));
867  }
868  }
870  dgdx_dot.getLinearOp() == Teuchos::null && dgdx_dot.isEmpty() == false,
871  std::logic_error,
872  "Error! Stokhos::SGModelEvaluator::evalModel: " <<
873  "Operator form of dg/dxdot(" << i << ") is required!");
874  }
875 
876  // dg/dx
877  if (outArgs.supports(OUT_ARG_DgDx, i).supports(DERIV_LINEAR_OP)) {
878  Derivative dgdx = outArgs.get_DgDx(i);
879  if (dgdx.getLinearOp() != Teuchos::null) {
881  Teuchos::rcp_dynamic_cast<Stokhos::SGOperator>(
882  dgdx.getLinearOp(), true);
884  op->getSGPolynomial();
885  if (me_outargs.supports(OUT_ARG_DgDx, ii).supports(DERIV_LINEAR_OP))
886  me_outargs.set_DgDx_sg(ii, sg_blocks);
887  else {
889  Teuchos::rcp_dynamic_cast<Stokhos::EpetraMultiVectorOperatorOrthogPoly>(sg_blocks, true);
891  sg_mv_op->multiVectorOrthogPoly();
892  if (me_outargs.supports(OUT_ARG_DgDx_sg, ii).supports(DERIV_MV_BY_COL))
893  me_outargs.set_DgDx_sg(ii, SGDerivative(dgdx_sg,
894  DERIV_MV_BY_COL));
895  else
896  me_outargs.set_DgDx_sg(ii, SGDerivative(dgdx_sg,
897  DERIV_TRANS_MV_BY_ROW));
898  }
899  }
901  dgdx.getLinearOp() == Teuchos::null && dgdx.isEmpty() == false,
902  std::logic_error,
903  "Error! Stokhos::SGModelEvaluator::evalModel: " <<
904  "Operator form of dg/dx(" << i << ") is required!");
905  }
906 
907  // dg/dp -- determinsitic p
908  for (int j=0; j<num_p; j++) {
909  if (!outArgs.supports(OUT_ARG_DgDp, i, j).none()) {
910  Derivative dgdp = outArgs.get_DgDp(i,j);
911  if (dgdp.getMultiVector() != Teuchos::null) {
913  if (dgdp.getMultiVectorOrientation() == DERIV_MV_BY_COL)
914  dgdp_sg =
916  sg_basis, overlapped_stoch_row_map,
917  me->get_g_map(ii), sg_g_map[i], sg_comm,
918  View, *(dgdp.getMultiVector())));
919  else if (dgdp.getMultiVectorOrientation() == DERIV_TRANS_MV_BY_ROW) {
921  Teuchos::rcp(&(dgdp.getMultiVector()->Map()),false);
922  dgdp_sg =
924  sg_basis, overlapped_stoch_row_map,
925  me->get_p_map(j), product_map, sg_comm,
926  View, *(dgdp.getMultiVector())));
927  }
928  me_outargs.set_DgDp_sg(
929  ii, j, SGDerivative(dgdp_sg, dgdp.getMultiVectorOrientation()));
930  }
931  else if (dgdp.getLinearOp() != Teuchos::null) {
933  Teuchos::rcp_dynamic_cast<Stokhos::EpetraOperatorOrthogPoly>(dgdp.getLinearOp(), true);
934  me_outargs.set_DgDp_sg(ii, j, SGDerivative(dgdp_sg));
935  }
936  }
937  }
938 
939  // dgdp -- stochastic p. Here we only support DERIV_LINEAR_OP
940  for (int j=0; j<num_p_sg; j++) {
941  if (!outArgs.supports(OUT_ARG_DgDp, i, j+num_p).none()) {
942  Derivative dgdp = outArgs.get_DgDp(i,j+num_p);
943  if (dgdp.getLinearOp() != Teuchos::null) {
945  Teuchos::rcp_dynamic_cast<Stokhos::MatrixFreeOperator>(dgdp.getLinearOp(), true);
947  dgdp_op->getSGPolynomial();
948  int jj = sg_p_index_map[j];
949  if (me_outargs.supports(OUT_ARG_DgDp_sg,ii,jj).supports(DERIV_LINEAR_OP))
950  me_outargs.set_DgDp_sg(ii, jj, SGDerivative(dgdp_op_sg));
951  else {
953  Teuchos::rcp_dynamic_cast<Stokhos::EpetraMultiVectorOperatorOrthogPoly>(dgdp_op_sg, true);
955  sg_mv_op->multiVectorOrthogPoly();
956  if (me_outargs.supports(OUT_ARG_DgDp_sg,ii,jj).supports(DERIV_MV_BY_COL))
957  me_outargs.set_DgDp_sg(
958  ii, jj, SGDerivative(dgdp_sg, DERIV_MV_BY_COL));
959  else
960  me_outargs.set_DgDp_sg(
961  ii, jj, SGDerivative(dgdp_sg, DERIV_TRANS_MV_BY_ROW));
962  }
963  }
965  dgdp.getLinearOp() == Teuchos::null && dgdp.isEmpty() == false,
966  std::logic_error,
967  "Error! Stokhos::SGModelEvaluator::evalModel: " <<
968  "Operator form of dg/dp(" << i << "," << j+num_p << ") is required!");
969  }
970  }
971 
972  }
973 
974  // Compute the functions
975  me->evalModel(me_inargs, me_outargs);
976 
977  // Copy block SG components for W
978  if ((W_out != Teuchos::null || (eval_W_with_f && f_out != Teuchos::null))
979  && !eval_prec) {
981  if (W_out != Teuchos::null)
982  W = W_out;
983  else
984  W = my_W;
986  Teuchos::rcp_dynamic_cast<Stokhos::SGOperator>(W, true);
987  W_sg->setupOperator(W_sg_blocks);
988 
989  if (WPrec_out != Teuchos::null) {
991  Teuchos::rcp_dynamic_cast<Stokhos::SGPreconditioner>(WPrec_out, true);
992  W_prec->setupPreconditioner(W_sg, *my_x);
993  }
994  }
995 
996  // f
997  if (f_out!=Teuchos::null){
998  if (!scaleOP)
999  for (int i=0; i<sg_basis->size(); i++)
1000  (*f_sg_blocks)[i].Scale(sg_basis->norm_squared(i));
1001  f_out->Export(*(f_sg_blocks->getBlockVector()), *sg_overlapped_f_exporter,
1002  Insert);
1003  }
1004 
1005  // df/dp -- deterministic p
1006  for (int i=0; i<num_p; i++) {
1007  if (!outArgs.supports(OUT_ARG_DfDp, i).none()) {
1008  Derivative dfdp = outArgs.get_DfDp(i);
1009  SGDerivative dfdp_sg = me_outargs.get_DfDp_sg(i);
1010  if (dfdp.getMultiVector() != Teuchos::null) {
1011  dfdp.getMultiVector()->Export(
1012  *(dfdp_sg.getMultiVector()->getBlockMultiVector()),
1013  *sg_overlapped_f_exporter, Insert);
1014  }
1015  // need to handle parallel df/dp operator case -- currently we have
1016  // no way to do communication of operators
1017  }
1018  }
1019 
1020  // df/dp -- stochastic p
1021  // No need to handle this case as it is handled by the operator
1022 
1023  // g, dg/dx, dg/dxdot, dg/dp -- these are all currently serial
1024 }
1025 
1026 void
1028  const Stokhos::EpetraVectorOrthogPoly& x_sg_in)
1029 {
1030  *sg_x_init = x_sg_in;
1031 }
1032 
1035 {
1036  return sg_x_init;
1037 }
1038 
1039 void
1041  int i, const Stokhos::EpetraVectorOrthogPoly& p_sg_in)
1042 {
1043  Teuchos::Array<int>::iterator it = std::find(sg_p_index_map.begin(),
1044  sg_p_index_map.end(),
1045  i);
1046  TEUCHOS_TEST_FOR_EXCEPTION(it == sg_p_index_map.end(), std::logic_error,
1047  "Error! Invalid p map index " << i);
1048  int ii = it - sg_p_index_map.begin();
1049  *sg_p_init[ii] = p_sg_in;
1050 }
1051 
1054 {
1055  Teuchos::Array<int>::const_iterator it = std::find(sg_p_index_map.begin(),
1056  sg_p_index_map.end(),
1057  l);
1058  TEUCHOS_TEST_FOR_EXCEPTION(it == sg_p_index_map.end(), std::logic_error,
1059  "Error! Invalid p map index " << l);
1060  int ll = it - sg_p_index_map.begin();
1061  return sg_p_init[ll];
1062 }
1063 
1066 {
1067  return sg_p_index_map;
1068 }
1069 
1072 {
1073  return sg_g_index_map;
1074 }
1075 
1078 {
1080  for (int i=0; i<num_g; i++)
1081  base_maps[i] = me->get_g_map(i);
1082  return base_maps;
1083  }
1084 
1087 {
1088  return overlapped_stoch_row_map;
1089 }
1090 
1093 {
1094  return sg_overlapped_x_map;
1095 }
1096 
1099 {
1100  return sg_overlapped_x_importer;
1101 }
1102 
1105  const Epetra_Vector* v) const
1106 {
1108  if (v == NULL)
1110  sg_basis, stoch_row_map, x_map, sg_x_map, sg_comm));
1111  else
1113  sg_basis, stoch_row_map, x_map, sg_x_map, sg_comm,
1114  CV, *v));
1115  return sg_x;
1116 }
1117 
1120  const Epetra_Vector* v) const
1121 {
1123  if (v == NULL)
1125  sg_basis, overlapped_stoch_row_map, x_map,
1126  sg_overlapped_x_map, sg_comm));
1127  else
1129  sg_basis, overlapped_stoch_row_map, x_map,
1130  sg_overlapped_x_map, sg_comm, CV, *v));
1131  return sg_x;
1132 }
1133 
1136  const Epetra_MultiVector* v) const
1137 {
1139  if (v == NULL)
1141  sg_basis, stoch_row_map, x_map, sg_x_map, sg_comm,
1142  num_vecs));
1143  else
1145  sg_basis, stoch_row_map, x_map, sg_x_map, sg_comm,
1146  CV, *v));
1147  return sg_x;
1148 }
1149 
1152  int num_vecs,
1153  Epetra_DataAccess CV,
1154  const Epetra_MultiVector* v) const
1155 {
1157  if (v == NULL)
1159  sg_basis, overlapped_stoch_row_map, x_map,
1160  sg_overlapped_x_map, sg_comm, num_vecs));
1161  else
1163  sg_basis, overlapped_stoch_row_map, x_map,
1164  sg_overlapped_x_map, sg_comm, CV, *v));
1165  return sg_x;
1166 }
1167 
1170  const Epetra_Vector* v) const
1171 {
1173  Teuchos::Array<int>::const_iterator it = std::find(sg_p_index_map.begin(),
1174  sg_p_index_map.end(),
1175  l);
1176  TEUCHOS_TEST_FOR_EXCEPTION(it == sg_p_index_map.end(), std::logic_error,
1177  "Error! Invalid p map index " << l);
1178  int ll = it - sg_p_index_map.begin();
1179  if (v == NULL)
1181  sg_basis, overlapped_stoch_p_map, me->get_p_map(l),
1182  sg_p_map[ll], sg_comm));
1183  else
1185  sg_basis, overlapped_stoch_p_map, me->get_p_map(l),
1186  sg_p_map[ll], sg_comm, CV, *v));
1187  return sg_p;
1188 }
1189 
1192  const Epetra_Vector* v) const
1193 {
1195  if (v == NULL)
1197  sg_basis, stoch_row_map, f_map, sg_f_map, sg_comm));
1198  else
1200  sg_basis, stoch_row_map, f_map, sg_f_map, sg_comm,
1201  CV, *v));
1202  return sg_f;
1203 }
1204 
1207  const Epetra_Vector* v) const
1208 {
1210  if (v == NULL)
1212  sg_basis, overlapped_stoch_row_map, f_map,
1213  sg_overlapped_f_map, sg_comm));
1214  else
1216  sg_basis, overlapped_stoch_row_map, f_map,
1217  sg_overlapped_f_map, sg_comm, CV, *v));
1218  return sg_f;
1219 }
1220 
1223  int num_vecs,
1224  Epetra_DataAccess CV,
1225  const Epetra_MultiVector* v) const
1226 {
1228  if (v == NULL)
1230  sg_basis, stoch_row_map, f_map, sg_f_map, sg_comm,
1231  num_vecs));
1232  else
1234  sg_basis, stoch_row_map, f_map, sg_f_map, sg_comm,
1235  CV, *v));
1236  return sg_f;
1237 }
1238 
1241  int num_vecs,
1242  Epetra_DataAccess CV,
1243  const Epetra_MultiVector* v) const
1244 {
1246  if (v == NULL)
1248  sg_basis, overlapped_stoch_row_map, f_map,
1249  sg_overlapped_f_map, sg_comm, num_vecs));
1250  else
1252  sg_basis, overlapped_stoch_row_map, f_map,
1253  sg_overlapped_f_map, sg_comm, CV, *v));
1254  return sg_f;
1255 }
1256 
1259  const Epetra_Vector* v) const
1260 {
1262  Teuchos::Array<int>::const_iterator it = std::find(sg_g_index_map.begin(),
1263  sg_g_index_map.end(),
1264  l);
1265  TEUCHOS_TEST_FOR_EXCEPTION(it == sg_g_index_map.end(), std::logic_error,
1266  "Error! Invalid g map index " << l);
1267  int ll = it - sg_g_index_map.begin();
1268  if (v == NULL)
1270  sg_basis, overlapped_stoch_row_map,
1271  me->get_g_map(l),
1272  sg_g_map[ll], sg_comm));
1273  else
1275  sg_basis, overlapped_stoch_row_map,
1276  me->get_g_map(l),
1277  sg_g_map[ll], sg_comm, CV, *v));
1278  return sg_g;
1279 }
1280 
1283  Epetra_DataAccess CV,
1284  const Epetra_MultiVector* v) const
1285 {
1287  Teuchos::Array<int>::const_iterator it = std::find(sg_g_index_map.begin(),
1288  sg_g_index_map.end(),
1289  l);
1290  TEUCHOS_TEST_FOR_EXCEPTION(it == sg_g_index_map.end(), std::logic_error,
1291  "Error! Invalid g map index " << l);
1292  int ll = it - sg_g_index_map.begin();
1293  if (v == NULL)
1295  sg_basis, overlapped_stoch_row_map,
1296  me->get_g_map(l),
1297  sg_g_map[ll], sg_comm, num_vecs));
1298  else
1300  sg_basis, overlapped_stoch_row_map,
1301  me->get_g_map(l),
1302  sg_g_map[ll], sg_comm, CV, *v));
1303  return sg_g;
1304 }
1305 
1308 {
1310  Teuchos::rcp(new EpetraExt::BlockVector(*x_map, *sg_overlapped_x_map));
1311  x_overlapped->Import(x, *sg_overlapped_x_importer, Insert);
1312  return x_overlapped;
1313 }
1314 
1317 {
1320  sg_basis, overlapped_stoch_row_map, x_map,
1321  sg_overlapped_x_map, sg_comm));
1322  sg_x_overlapped->getBlockVector()->Import(x, *sg_overlapped_x_importer,
1323  Insert);
1324  return sg_x_overlapped;
1325 }
1326 
1329 {
1331  Teuchos::rcp(new EpetraExt::BlockVector(*x_map, *sg_x_map));
1332  x->Export(x_overlapped, *sg_overlapped_x_importer, Insert);
1333  return x;
1334 }
1335 
1338 {
1340  Teuchos::rcp(new EpetraExt::BlockVector(*f_map, *sg_overlapped_f_map));
1341  f_overlapped->Import(f, *sg_overlapped_f_exporter, Insert);
1342  return f_overlapped;
1343 }
1344 
1347 {
1349  Teuchos::rcp(new EpetraExt::BlockVector(*f_map, *sg_f_map));
1350  f->Export(f_overlapped, *sg_overlapped_f_exporter, Insert);
1351  return f;
1352 }
Teuchos::RCP< const Stokhos::EpetraSparse3Tensor > epetraCijk
Epetra Cijk.
Teuchos::RCP< Stokhos::EpetraVectorOrthogPoly > create_x_sg(Epetra_DataAccess CV=Copy, const Epetra_Vector *v=NULL) const
Create vector orthog poly using x map and owned sg map.
Teuchos::RCP< Stokhos::EpetraVectorOrthogPoly > create_p_sg(int l, Epetra_DataAccess CV=Copy, const Epetra_Vector *v=NULL) const
Create vector orthog poly using p map.
Teuchos::RCP< Stokhos::EpetraOperatorOrthogPoly > W_sg_blocks
W stochastic Galerkin components.
Teuchos::RCP< const Epetra_Map > get_x_map() const
Return solution vector map.
Teuchos::RCP< Stokhos::EpetraMultiVectorOrthogPoly > create_f_mv_sg_overlap(int num_vecs, Epetra_DataAccess CV=Copy, const Epetra_MultiVector *v=NULL) const
Create multi-vector orthog poly using f map and overlap sg map.
An Epetra operator representing the block stochastic Galerkin operator.
Teuchos::RCP< const Epetra_Map > sg_overlapped_f_map
Block SG overlapped residual map.
void setCoeffPtr(ordinal_type i, const Teuchos::RCP< coeff_type > &c)
Set coefficient i to c.
virtual void setupOperator(const Teuchos::RCP< Stokhos::EpetraOperatorOrthogPoly > &poly)=0
Setup operator.
Teuchos::RCP< Stokhos::EpetraVectorOrthogPoly > sg_x_init
SG initial x.
A container class storing an orthogonal polynomial whose coefficients are vectors, operators, or in general any type that would have an expensive copy constructor.
Teuchos::Array< int > sg_g_index_map
Index map between block-g and g_sg maps.
unsigned int num_sg_blocks
Number of stochastic blocks.
Teuchos::RCP< const Epetra_Map > get_g_map(int l) const
Return response map.
bool myGID(int i) const
Return whether global index i resides on this processor.
Teuchos::RCP< const Epetra_Map > x_map
Underlying unknown map.
Teuchos::RCP< const Epetra_Map > get_f_map() const
Return residual vector map.
Teuchos::RCP< Stokhos::EpetraVectorOrthogPoly > x_sg_blocks
x stochastic Galerkin components
Teuchos::RCP< Stokhos::EpetraVectorOrthogPoly > f_sg_blocks
f stochastic Galerkin components
int num_p_sg
Number of stochastic parameter vectors.
T & get(ParameterList &l, const std::string &name)
bool eval_W_with_f
Whether to always evaluate W with f.
void evalModel(const InArgs &inArgs, const OutArgs &outArgs) const
Evaluate model on InArgs.
Teuchos::RCP< EpetraExt::BlockVector > getBlockVector()
Get block vector.
Teuchos::RCP< const Epetra_BlockMap > get_overlap_stochastic_map() const
Return overlap stochastic map.
Teuchos::RCP< const Epetra_Map > sg_overlapped_x_map
Block SG overlapped unknown map.
SGModelEvaluator(const Teuchos::RCP< EpetraExt::ModelEvaluator > &me, const Teuchos::RCP< const Stokhos::OrthogPolyBasis< int, double > > &sg_basis, const Teuchos::RCP< const Stokhos::Quadrature< int, double > > &sg_quad, const Teuchos::RCP< Stokhos::OrthogPolyExpansion< int, double > > &sg_exp, const Teuchos::RCP< const Stokhos::ParallelData > &sg_parallel_data, const Teuchos::RCP< Teuchos::ParameterList > &params, bool scaleOP=true)
#define TEUCHOS_TEST_FOR_EXCEPTION(throw_exception_test, Exception, msg)
Teuchos::Array< Teuchos::RCP< const Epetra_Map > > sg_g_map
Block SG response map.
Teuchos::RCP< const Stokhos::EpetraVectorOrthogPoly > get_p_sg_init(int l) const
Return initial SG parameters.
Teuchos::RCP< Epetra_Operator > create_DgDp_op(int j, int i) const
Create SG operator representing dg/dp.
unsigned int num_p_blocks
Number of p stochastic blocks (may be smaller than num_sg_blocks)
RCP< ParameterList > sublist(const RCP< ParameterList > &paramList, const std::string &name, bool mustAlreadyExist=false, const std::string &docString="")
T * get() const
ParameterList & set(std::string const &name, T &&value, std::string const &docString="", RCP< const ParameterEntryValidator > const &validator=null)
Teuchos::RCP< const Stokhos::ParallelData > sg_parallel_data
Parallel SG data.
Teuchos::RCP< const Epetra_BlockMap > overlapped_stoch_p_map
Overlapped map for p stochastic blocks (local map)
virtual ordinal_type dimension() const =0
Return dimension of basis.
Teuchos::Array< Teuchos::RCP< Stokhos::EpetraVectorOrthogPoly > > sg_p_init
SG initial p.
Teuchos::RCP< const Stokhos::EpetraSparse3Tensor > serialCijk
Serial Epetra Cijk for dgdx*.
Teuchos::Array< int > get_p_sg_map_indices() const
Get indices of SG parameters.
Teuchos::RCP< const Epetra_Map > sg_f_map
Block SG residual map.
OutArgs createOutArgs() const
Create OutArgs.
Abstract base class for orthogonal polynomial-based expansions.
Teuchos::RCP< EpetraExt::BlockVector > import_solution(const Epetra_Vector &x) const
Import parallel solution vector.
int NumMyElements() const
Teuchos::RCP< Stokhos::EpetraVectorOrthogPoly > create_g_sg(int l, Epetra_DataAccess CV=Copy, const Epetra_Vector *v=NULL) const
Create vector orthog poly using g map.
Teuchos::RCP< Stokhos::EpetraMultiVectorOrthogPoly > create_f_mv_sg(int num_vecs, Epetra_DataAccess CV=Copy, const Epetra_MultiVector *v=NULL) const
Create multi-vector orthog poly using f map and owned sg map.
int num_g_sg
Number of stochastic response vectors.
Teuchos::RCP< const Epetra_BlockMap > stoch_row_map
Map for stochastic blocks.
Teuchos::RCP< const Epetra_BlockMap > get_x_sg_overlap_map() const
Return x sg overlap map.
TEUCHOS_DEPRECATED RCP< T > rcp(T *p, Dealloc_T dealloc, bool owns_mem)
An abstract class to represent a generic stochastic Galerkin operator as an Epetra_Operator.
Teuchos::RCP< Stokhos::EpetraVectorOrthogPoly > x_dot_sg_blocks
x_dot stochastic Galerkin components
Teuchos::RCP< const Epetra_Vector > get_x_init() const
Return initial solution.
Teuchos::RCP< Stokhos::SGPreconditionerFactory > sg_prec_factory
Preconditioner factory.
A container class storing an orthogonal polynomial whose coefficients are vectors, operators, or in general any type that would have an expensive copy constructor.
virtual Teuchos::RCP< Stokhos::SGOperator > build(const Teuchos::RCP< const EpetraExt::MultiComm > &sg_comm, const Teuchos::RCP< const Stokhos::OrthogPolyBasis< int, double > > &sg_basis, const Teuchos::RCP< const Stokhos::EpetraSparse3Tensor > &epetraCijk, const Teuchos::RCP< const Epetra_Map > &domain_base_map, const Teuchos::RCP< const Epetra_Map > &range_base_map, const Teuchos::RCP< const Epetra_Map > &domain_sg_map, const Teuchos::RCP< const Epetra_Map > &range_sg_map)
Build preconditioner operator.
virtual Teuchos::RCP< Stokhos::EpetraOperatorOrthogPoly > getSGPolynomial()=0
Get SG polynomial.
Teuchos::RCP< const Stokhos::OrthogPolyBasis< int, double > > sg_basis
Stochastic Galerkin basis.
Teuchos::Array< Teuchos::RCP< const Epetra_Map > > sg_p_map
Block SG parameter map.
A container class storing an orthogonal polynomial whose coefficients are vectors, operators, or in general any type that would have an expensive copy constructor.
bool supports_x
Whether we support x (and thus f and W)
void set_p_sg_init(int i, const Stokhos::EpetraVectorOrthogPoly &p_sg_in)
Set initial parameter polynomial.
Teuchos::RCP< const Stokhos::EpetraVectorOrthogPoly > get_x_sg_init() const
Return initial SG x.
Teuchos::RCP< EpetraMultiVectorOrthogPoly > multiVectorOrthogPoly() const
Get multi vector orthog poly.
Teuchos::Array< int > get_g_sg_map_indices() const
Get indices of SG responses.
Factory for generating stochastic Galerkin preconditioners.
Teuchos::RCP< Epetra_Export > sg_overlapped_f_exporter
Exporter from SG-overlapped to SG maps.
Teuchos::RCP< Stokhos::EpetraVectorOrthogPoly > create_x_sg_overlap(Epetra_DataAccess CV=Copy, const Epetra_Vector *v=NULL) const
Create vector orthog poly using x map and overlap sg map.
virtual Teuchos::RCP< Stokhos::EpetraOperatorOrthogPoly > getSGPolynomial()
Get SG polynomial.
Teuchos::RCP< Stokhos::EpetraMultiVectorOrthogPoly > create_x_mv_sg(int num_vecs, Epetra_DataAccess CV=Copy, const Epetra_MultiVector *v=NULL) const
Create vector orthog poly using x map and owned sg map.
Teuchos::Array< Teuchos::RCP< Teuchos::Array< std::string > > > sg_p_names
SG coefficient parameter names.
Teuchos::RCP< EpetraExt::ModelEvaluator > me
Underlying model evaluator.
Teuchos::RCP< const Epetra_Map > get_p_map(int l) const
Return parameter vector map.
Teuchos::RCP< Stokhos::EpetraVectorOrthogPoly > create_f_sg_overlap(Epetra_DataAccess CV=Copy, const Epetra_Vector *v=NULL) const
Create vector orthog poly using f map and overlap sg map.
iterator end()
std::vector< T >::const_iterator const_iterator
virtual void setupPreconditioner(const Teuchos::RCP< Stokhos::SGOperator > &sg_op, const Epetra_Vector &x)=0
Setup preconditioner.
Teuchos::RCP< Stokhos::EpetraMultiVectorOrthogPoly > create_g_mv_sg(int l, int num_vecs, Epetra_DataAccess CV=Copy, const Epetra_MultiVector *v=NULL) const
Create multi-vector orthog poly using g map.
unsigned int num_W_blocks
Number of W stochastic blocks (may be smaller than num_sg_blocks)
InArgs createInArgs() const
Create InArgs.
void set_x_sg_init(const Stokhos::EpetraVectorOrthogPoly &x_sg_in)
Set initial solution polynomial.
Teuchos::RCP< Stokhos::EpetraVectorOrthogPoly > import_solution_poly(const Epetra_Vector &x) const
Import parallel solution vector.
Teuchos::RCP< EpetraExt::BlockVector > import_residual(const Epetra_Vector &f) const
Import parallel residual vector.
Teuchos::RCP< const Epetra_Import > get_x_sg_importer() const
Return x sg importer.
void push_back(const value_type &x)
Teuchos::RCP< const Epetra_BlockMap > overlapped_stoch_row_map
Overlapped map for stochastic blocks (local map)
ScalarType f(const Teuchos::Array< ScalarType > &x, double a, double b)
int num_g
Number of response vectors of underlying model evaluator.
Teuchos::RCP< Epetra_Operator > create_DgDx_dot_op(int j) const
Create SG operator representing dg/dxdot.
int num_p
Number of parameter vectors of underlying model evaluator.
Teuchos::RCP< Epetra_Import > sg_overlapped_x_importer
Importer from SG to SG-overlapped maps.
Teuchos::RCP< Epetra_Operator > create_DgDx_op(int j) const
Create SG operator representing dg/dx.
size_type size() const
Factory for generating stochastic Galerkin preconditioners.
Teuchos::RCP< Epetra_Vector > my_x
x pointer for evaluating preconditioner
Teuchos::RCP< Teuchos::ParameterList > params
Algorithmic parameters.
Teuchos::RCP< const Teuchos::Array< std::string > > get_p_names(int l) const
Return array of parameter names.
Teuchos::RCP< EpetraExt::BlockVector > export_residual(const Epetra_Vector &f_overlapped) const
Export parallel residual vector.
Teuchos::RCP< EpetraExt::BlockVector > export_solution(const Epetra_Vector &x_overlapped) const
Export parallel solution vector.
Teuchos::RCP< Epetra_Operator > create_W() const
Create W = alpha*M + beta*J matrix.
const Epetra_CrsGraph & Graph() const
An abstract class to represent a generic stochastic Galerkin preconditioner as an Epetra_Operator...
Teuchos::Array< int > sg_p_index_map
Index map between block-p and p_sg maps.
Epetra_DataAccess
Teuchos::Array< Teuchos::RCP< const Epetra_Map > > get_g_sg_base_maps() const
Get base maps of SG responses.
Teuchos::RCP< const Epetra_Map > f_map
Underlying residual map.
std::vector< T >::iterator iterator
iterator begin()
Teuchos::RCP< Stokhos::EpetraMultiVectorOrthogPoly > create_x_mv_sg_overlap(int num_vecs, Epetra_DataAccess CV=Copy, const Epetra_MultiVector *v=NULL) const
Create vector orthog poly using x map and overlap sg map.
Teuchos::RCP< Stokhos::EpetraVectorOrthogPoly > create_f_sg(Epetra_DataAccess CV=Copy, const Epetra_Vector *v=NULL) const
Create vector orthog poly using f map and owned sg map.
Teuchos::RCP< const Epetra_Vector > get_p_init(int l) const
Return initial parameters.
Teuchos::RCP< EpetraExt::ModelEvaluator::Preconditioner > create_WPrec() const
Create preconditioner operator.
Teuchos::RCP< Epetra_Operator > create_DfDp_op(int i) const
Create SG operator representing df/dp.
ScalarType g(const Teuchos::Array< ScalarType > &x, const ScalarType &y)
Teuchos::RCP< const Epetra_Map > sg_x_map
Block SG unknown map.
Teuchos::RCP< const EpetraExt::MultiComm > sg_comm
Parallel SG communicator.
A container class storing an orthogonal polynomial whose coefficients are vectors, operators, or in general any type that would have an expensive copy constructor.