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_Adaptive.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 
27  const Teuchos::RCP<const Stokhos::Quadrature<int,double> >& sg_quad_,
29  const Teuchos::RCP<const Stokhos::ParallelData>& sg_parallel_data_,
30  bool onlyUseLinear_,int kExpOrder_,
32  : me(me_),
33  sg_basis(am->getMasterStochasticBasis()),
34  sg_row_dof_basis(am->getRowStochasticBasis()),
35  sg_quad(sg_quad_),
36  sg_exp(sg_exp_),
37  params(params_),
38  num_sg_blocks(sg_basis->size()),
39  num_W_blocks(sg_basis->size()),
40  num_p_blocks(sg_basis->size()),
41  supports_x(false),
42  x_map(me->get_x_map()),
43  f_map(me->get_f_map()),
44  sg_parallel_data(sg_parallel_data_),
45  sg_comm(sg_parallel_data->getMultiComm()),
46  epetraCijk(sg_parallel_data->getEpetraCijk()),
47  serialCijk(),
48  Cijk(epetraCijk->getParallelCijk()),
49  num_p(0),
50  num_p_sg(0),
51  sg_p_index_map(),
52  sg_p_map(),
53  sg_p_names(),
54  num_g(0),
55  num_g_sg(0),
56  sg_g_index_map(),
57  sg_g_map(),
58  x_dot_sg_blocks(),
59  x_sg_blocks(),
60  f_sg_blocks(),
61  W_sg_blocks(),
62  dgdx_dot_sg_blocks(),
63  dgdx_sg_blocks(),
64  sg_x_init(),
65  sg_p_init(),
66  eval_W_with_f(false),
67  kExpOrder(kExpOrder_),
68  onlyUseLinear(onlyUseLinear_),
69  scaleOP(am->isScaled()),
70  adaptMngr(am)
71 {
72  if (x_map != Teuchos::null)
73  supports_x = true;
74 
76  Teuchos::rcp(new Epetra_LocalMap(static_cast<int>(num_sg_blocks), 0, *(sg_parallel_data->getStochasticComm())));
78 
80  stoch_row_map = epetraCijk->getStochasticRowMap();
81 
82  if (supports_x) {
83 
84  // Create block SG x and f maps
85  sg_x_map =
86  Teuchos::rcp(EpetraExt::BlockUtility::GenerateBlockMap(
87  *x_map, *stoch_row_map, *sg_comm));
89  Teuchos::rcp(EpetraExt::BlockUtility::GenerateBlockMap(
91 
92  sg_f_map =
93  Teuchos::rcp(EpetraExt::BlockUtility::GenerateBlockMap(
94  *f_map, *stoch_row_map, *sg_comm));
96  Teuchos::rcp(EpetraExt::BlockUtility::GenerateBlockMap(
98 
99  // Create interlace SG x and f maps
100  adapted_x_map = adaptMngr->getAdaptedMap();
102 
103  adapted_f_map = adapted_x_map; // these are the same! No overlap possible in refinement!
105 
106  // Create importer/exporter from/to overlapped distribution
111 
112  // now we create the underlying Epetra block vectors
113  // that will be used by the model evaluator to construct
114  // the solution of the deterministic problem.
116 
117  // Create vector blocks
121 
122  // Create default sg_x_init
124  if (sg_x_init->myGID(0))
125  (*sg_x_init)[0] = *(me->get_x_init());
126 
127  // Preconditioner needs an x: This is adapted
129 
130  // setup storage for W, these are blocked in Stokhos
131  // format
133 
134  // Determine W expansion type
135  std::string W_expansion_type =
136  params->get("Jacobian Expansion Type", "Full");
137  if (W_expansion_type == "Linear")
139  else
141 
142  Teuchos::RCP<Epetra_BlockMap> W_overlap_map =
144  static_cast<int>(num_W_blocks), 0,
145  *(sg_parallel_data->getStochasticComm())));
146  W_sg_blocks =
148  sg_basis, W_overlap_map, x_map, f_map, adapted_f_map,
149  sg_comm));
150  for (unsigned int i=0; i<num_W_blocks; i++)
151  W_sg_blocks->setCoeffPtr(i, me->create_W()); // allocate a bunch of matrices
152 
153  eval_W_with_f = params->get("Evaluate W with F", false);
154  }
155 
156  // Parameters -- The idea here is to add new parameter vectors
157  // for the stochastic Galerkin components of the parameters
158 
159  InArgs me_inargs = me->createInArgs();
160  OutArgs me_outargs = me->createOutArgs();
161  num_p = me_inargs.Np();
162 
163  // Get the p_sg's supported and build index map
164  for (int i=0; i<num_p; i++) {
165  if (me_inargs.supports(IN_ARG_p_sg, i))
167  }
169 
170  sg_p_map.resize(num_p_sg);
171  sg_p_names.resize(num_p_sg);
172  sg_p_init.resize(num_p_sg);
173 
174  // Determine parameter expansion type
175  std::string p_expansion_type =
176  params->get("Parameter Expansion Type", "Full");
177  if (p_expansion_type == "Linear")
179  else
181 
182  // Create parameter maps, names, and initial values
185  static_cast<int>(num_p_blocks), 0,
186  *(sg_parallel_data->getStochasticComm())));
187  for (int i=0; i<num_p_sg; i++) {
188  Teuchos::RCP<const Epetra_Map> p_map = me->get_p_map(sg_p_index_map[i]);
189  sg_p_map[i] =
190  Teuchos::rcp(EpetraExt::BlockUtility::GenerateBlockMap(
191  *p_map, *overlapped_stoch_p_map, *sg_comm));
192 
194  me->get_p_names(sg_p_index_map[i]);
195  if (p_names != Teuchos::null) {
196  sg_p_names[i] =
198  for (int j=0; j<p_names->size(); j++) {
199  std::stringstream ss;
200  ss << (*p_names)[j] << " -- SG Coefficient " << i;
201  (*sg_p_names[i])[j] = ss.str();
202  }
203  }
204 
205  // Create default sg_p_init
207  sg_p_init[i]->init(0.0);
208  }
209 
210  // Responses -- The idea here is to add new parameter vectors
211  // for the stochastic Galerkin components of the respones
212  num_g = me_outargs.Ng();
213 
214  // Get the g_sg's supported and build index map
215  for (int i=0; i<num_g; i++) {
216  if (me_outargs.supports(OUT_ARG_g_sg, i))
218  }
220 
221  sg_g_map.resize(num_g_sg);
223  dgdx_sg_blocks.resize(num_g_sg);
224 
225  // Create response maps
226  for (int i=0; i<num_g_sg; i++) {
227  Teuchos::RCP<const Epetra_Map> g_map = me->get_g_map(sg_g_index_map[i]);
228  sg_g_map[i] =
229  Teuchos::rcp(EpetraExt::BlockUtility::GenerateBlockMap(
230  *g_map, *overlapped_stoch_row_map, *sg_comm));
231 
232  // Create dg/dxdot, dg/dx SG blocks
233  if (supports_x) {
234  dgdx_dot_sg_blocks[i] =
237  dgdx_sg_blocks[i] =
240  }
241  }
242 
243  // We don't support parallel for dgdx yet, so build a new EpetraCijk
244  if (supports_x) {
245  serialCijk =
247  epetraCijk->getCijk(),
248  sg_comm,
250  }
251 
252 }
253 
256  const Teuchos::RCP<const Stokhos::OrthogPolyBasis<int,double> >& sg_basis_,
257  const std::vector<Teuchos::RCP<const Stokhos::ProductBasis<int,double> > > & sg_row_dof_basis_,
258  const Teuchos::RCP<const Stokhos::Quadrature<int,double> >& sg_quad_,
260  const Teuchos::RCP<const Stokhos::ParallelData>& sg_parallel_data_,
261  bool onlyUseLinear_,int kExpOrder_,
263  bool scaleOP_)
264  : me(me_),
265  sg_basis(sg_basis_),
266  sg_row_dof_basis(sg_row_dof_basis_),
267  sg_quad(sg_quad_),
268  sg_exp(sg_exp_),
269  params(params_),
270  num_sg_blocks(sg_basis->size()),
271  num_W_blocks(sg_basis->size()),
272  num_p_blocks(sg_basis->size()),
273  supports_x(false),
274  x_map(me->get_x_map()),
275  f_map(me->get_f_map()),
276  sg_parallel_data(sg_parallel_data_),
277  sg_comm(sg_parallel_data->getMultiComm()),
278  epetraCijk(sg_parallel_data->getEpetraCijk()),
279  serialCijk(),
280  Cijk(epetraCijk->getParallelCijk()),
281  num_p(0),
282  num_p_sg(0),
283  sg_p_index_map(),
284  sg_p_map(),
285  sg_p_names(),
286  num_g(0),
287  num_g_sg(0),
288  sg_g_index_map(),
289  sg_g_map(),
290  x_dot_sg_blocks(),
291  x_sg_blocks(),
292  f_sg_blocks(),
293  W_sg_blocks(),
294  dgdx_dot_sg_blocks(),
295  dgdx_sg_blocks(),
296  sg_x_init(),
297  sg_p_init(),
298  eval_W_with_f(false),
299  kExpOrder(kExpOrder_),
300  onlyUseLinear(onlyUseLinear_),
301  scaleOP(scaleOP_)
302 {
303  if (x_map != Teuchos::null)
304  supports_x = true;
305 
307  sg_row_dof_basis,*sg_parallel_data->getStochasticComm(),scaleOP));
308 
310  Teuchos::rcp(new Epetra_LocalMap(static_cast<int>(num_sg_blocks), 0, *(sg_parallel_data->getStochasticComm())));
312 
313  if (epetraCijk != Teuchos::null)
314  stoch_row_map = epetraCijk->getStochasticRowMap();
315 
316  if (supports_x) {
317 
318  // Create block SG x and f maps
319  sg_x_map =
320  Teuchos::rcp(EpetraExt::BlockUtility::GenerateBlockMap(
321  *x_map, *stoch_row_map, *sg_comm));
323  Teuchos::rcp(EpetraExt::BlockUtility::GenerateBlockMap(
325 
326  sg_f_map =
327  Teuchos::rcp(EpetraExt::BlockUtility::GenerateBlockMap(
328  *f_map, *stoch_row_map, *sg_comm));
330  Teuchos::rcp(EpetraExt::BlockUtility::GenerateBlockMap(
332 
333  // Create interlace SG x and f maps
334  adapted_x_map = adaptMngr->getAdaptedMap();
336 
337  adapted_f_map = adapted_x_map; // these are the same! No overlap possible in refinement!
339 
340  // Create importer/exporter from/to overlapped distribution
345 
346  // now we create the underlying Epetra block vectors
347  // that will be used by the model evaluator to construct
348  // the solution of the deterministic problem.
350 
351  // Create vector blocks
355 
356  // Create default sg_x_init
358  if (sg_x_init->myGID(0))
359  (*sg_x_init)[0] = *(me->get_x_init());
360 
361  // Preconditioner needs an x: This is adapted
363 
364  // setup storage for W, these are blocked in Stokhos
365  // format
367 
368  // Determine W expansion type
369  std::string W_expansion_type =
370  params->get("Jacobian Expansion Type", "Full");
371  if (W_expansion_type == "Linear")
373  else
375 
376  Teuchos::RCP<Epetra_BlockMap> W_overlap_map =
378  static_cast<int>(num_W_blocks), 0,
379  *(sg_parallel_data->getStochasticComm())));
380  W_sg_blocks =
382  sg_basis, W_overlap_map, x_map, f_map, adapted_f_map,
383  sg_comm));
384  for (unsigned int i=0; i<num_W_blocks; i++)
385  W_sg_blocks->setCoeffPtr(i, me->create_W()); // allocate a bunch of matrices
386 
387  eval_W_with_f = params->get("Evaluate W with F", false);
388  }
389 
390  // Parameters -- The idea here is to add new parameter vectors
391  // for the stochastic Galerkin components of the parameters
392 
393  InArgs me_inargs = me->createInArgs();
394  OutArgs me_outargs = me->createOutArgs();
395  num_p = me_inargs.Np();
396 
397  // Get the p_sg's supported and build index map
398  for (int i=0; i<num_p; i++) {
399  if (me_inargs.supports(IN_ARG_p_sg, i))
401  }
403 
404  sg_p_map.resize(num_p_sg);
405  sg_p_names.resize(num_p_sg);
406  sg_p_init.resize(num_p_sg);
407 
408  // Determine parameter expansion type
409  std::string p_expansion_type =
410  params->get("Parameter Expansion Type", "Full");
411  if (p_expansion_type == "Linear")
413  else
415 
416  // Create parameter maps, names, and initial values
419  static_cast<int>(num_p_blocks), 0,
420  *(sg_parallel_data->getStochasticComm())));
421  for (int i=0; i<num_p_sg; i++) {
422  Teuchos::RCP<const Epetra_Map> p_map = me->get_p_map(sg_p_index_map[i]);
423  sg_p_map[i] =
424  Teuchos::rcp(EpetraExt::BlockUtility::GenerateBlockMap(
425  *p_map, *overlapped_stoch_p_map, *sg_comm));
426 
428  me->get_p_names(sg_p_index_map[i]);
429  if (p_names != Teuchos::null) {
430  sg_p_names[i] =
432  for (int j=0; j<p_names->size(); j++) {
433  std::stringstream ss;
434  ss << (*p_names)[j] << " -- SG Coefficient " << i;
435  (*sg_p_names[i])[j] = ss.str();
436  }
437  }
438 
439  // Create default sg_p_init
441  sg_p_init[i]->init(0.0);
442  }
443 
444  // Responses -- The idea here is to add new parameter vectors
445  // for the stochastic Galerkin components of the respones
446  num_g = me_outargs.Ng();
447 
448  // Get the g_sg's supported and build index map
449  for (int i=0; i<num_g; i++) {
450  if (me_outargs.supports(OUT_ARG_g_sg, i))
452  }
454 
455  sg_g_map.resize(num_g_sg);
457  dgdx_sg_blocks.resize(num_g_sg);
458 
459  // Create response maps
460  for (int i=0; i<num_g_sg; i++) {
461  Teuchos::RCP<const Epetra_Map> g_map = me->get_g_map(sg_g_index_map[i]);
462  sg_g_map[i] =
463  Teuchos::rcp(EpetraExt::BlockUtility::GenerateBlockMap(
464  *g_map, *overlapped_stoch_row_map, *sg_comm));
465 
466  // Create dg/dxdot, dg/dx SG blocks
467  if (supports_x) {
468  dgdx_dot_sg_blocks[i] =
471  dgdx_sg_blocks[i] =
474  }
475  }
476 
477  // We don't support parallel for dgdx yet, so build a new EpetraCijk
478  if (supports_x) {
479  serialCijk =
481  epetraCijk->getCijk(),
482  sg_comm,
484  }
485 
486 }
487 
488 // Overridden from EpetraExt::ModelEvaluator
489 
492 {
493  return adapted_x_map;
494 }
495 
498 {
499  return adapted_f_map;
500 }
501 
504 {
505  TEUCHOS_TEST_FOR_EXCEPTION(l < 0 || l >= num_p + num_p_sg, std::logic_error,
506  "Error! Invalid p map index " << l);
507  if (l < num_p)
508  return me->get_p_map(l);
509  else
510  return sg_p_map[l-num_p];
511 
512  return Teuchos::null;
513 }
514 
517 {
518  TEUCHOS_TEST_FOR_EXCEPTION(l < 0 || l >= num_g_sg, std::logic_error,
519  "Error! Invalid g map index " << l);
520  return sg_g_map[l];
521 }
522 
525 {
526  TEUCHOS_TEST_FOR_EXCEPTION(l < 0 || l >= num_p + num_p_sg, std::logic_error,
527  "Error! Invalid p map index " << l);
528  if (l < num_p)
529  return me->get_p_names(l);
530  else
531  return sg_p_names[l-num_p];
532 
533  return Teuchos::null;
534 }
535 
538 {
539  // get stochastic galerking initial condition and write it out to x initial condition
540  Teuchos::RCP<Epetra_Vector> x_init = Teuchos::rcp(new Epetra_Vector(*get_x_map()));
541  adaptMngr->copyToAdaptiveVector(*get_x_sg_init(),*x_init);
542 
543  return x_init;
544 }
545 
548 {
549  TEUCHOS_TEST_FOR_EXCEPTION(l < 0 || l >= num_p + num_p_sg, std::logic_error,
550  "Error! Invalid p map index " << l);
551  if (l < num_p)
552  return me->get_p_init(l);
553  else
554  return sg_p_init[l-num_p]->getBlockVector();
555 
556  return Teuchos::null;
557 }
558 
561 {
562  if (supports_x) {
564  Teuchos::rcp(&(params->sublist("SG Operator")), false);
565 
567  = Teuchos::rcp_dynamic_cast<Epetra_CrsMatrix>(me->create_W(), true);
568  const Epetra_CrsGraph & W_graph = W_crs->Graph();
569 
570  adaptMngr->setupWithGraph(W_graph,onlyUseLinear,kExpOrder);
571  my_W = adaptMngr->buildMatrixFromGraph();
572  adaptMngr->setupOperator(*my_W,*Cijk,*W_sg_blocks);
573 
574  return my_W;
575  }
576 
577  return Teuchos::null;
578 }
579 
580 EpetraExt::ModelEvaluator::InArgs
582 {
583  InArgsSetup inArgs;
584  InArgs me_inargs = me->createInArgs();
585 
586  inArgs.setModelEvalDescription(this->description());
587  inArgs.set_Np(num_p + num_p_sg);
588  inArgs.setSupports(IN_ARG_x_dot, me_inargs.supports(IN_ARG_x_dot_sg));
589  inArgs.setSupports(IN_ARG_x, me_inargs.supports(IN_ARG_x_sg));
590  inArgs.setSupports(IN_ARG_t, me_inargs.supports(IN_ARG_t));
591  inArgs.setSupports(IN_ARG_alpha, me_inargs.supports(IN_ARG_alpha));
592  inArgs.setSupports(IN_ARG_beta, me_inargs.supports(IN_ARG_beta));
593  inArgs.setSupports(IN_ARG_sg_basis, me_inargs.supports(IN_ARG_sg_basis));
594  inArgs.setSupports(IN_ARG_sg_quadrature,
595  me_inargs.supports(IN_ARG_sg_quadrature));
596  inArgs.setSupports(IN_ARG_sg_expansion,
597  me_inargs.supports(IN_ARG_sg_expansion));
598 
599  return inArgs;
600 }
601 
602 EpetraExt::ModelEvaluator::OutArgs
604 {
605  OutArgsSetup outArgs;
606  OutArgs me_outargs = me->createOutArgs();
607 
608  outArgs.setModelEvalDescription(this->description());
609  outArgs.set_Np_Ng(num_p+num_p_sg, num_g_sg);
610 
611  outArgs.setSupports(OUT_ARG_f, me_outargs.supports(OUT_ARG_f_sg));
612  outArgs.setSupports(OUT_ARG_W, me_outargs.supports(OUT_ARG_W_sg));
613  outArgs.setSupports(OUT_ARG_WPrec, false);
614  for (int j=0; j<num_p; j++)
615  outArgs.setSupports(OUT_ARG_DfDp, j,
616  me_outargs.supports(OUT_ARG_DfDp_sg, j));
617 
618  for (int i=0; i<num_g_sg; i++) {
619  int ii = sg_g_index_map[i];
620 // if (!me_outargs.supports(OUT_ARG_DgDx_dot_sg, ii).none())
621 // outArgs.setSupports(OUT_ARG_DgDx_dot, i, DERIV_LINEAR_OP);
622 // if (!me_outargs.supports(OUT_ARG_DgDx_sg, ii).none())
623 // outArgs.setSupports(OUT_ARG_DgDx, i, DERIV_LINEAR_OP);
624  for (int j=0; j<num_p; j++)
625  outArgs.setSupports(OUT_ARG_DgDp, i, j,
626  me_outargs.supports(OUT_ARG_DgDp_sg, ii, j));
627  }
628 
629  // We do not support derivatives w.r.t. the new SG parameters, so their
630  // support defaults to none.
631 
632  return outArgs;
633 }
634 
635 void
637  const OutArgs& outArgs) const
638 {
639  // Get the input arguments
641  if (inArgs.supports(IN_ARG_x)) {
642  x = inArgs.get_x();
643  if (x != Teuchos::null)
644  *my_x = *x;
645  }
647  if (inArgs.supports(IN_ARG_x_dot))
648  x_dot = inArgs.get_x_dot();
649 
650  // Get the output arguments
651  EpetraExt::ModelEvaluator::Evaluation<Epetra_Vector> f_out;
652  if (outArgs.supports(OUT_ARG_f))
653  f_out = outArgs.get_f();
655  if (outArgs.supports(OUT_ARG_W))
656  W_out = outArgs.get_W();
657 
658  // Create underlying inargs
659  InArgs me_inargs = me->createInArgs();
660  if (x != Teuchos::null) {
661  // copy to a poly orthog vector
662  adaptMngr->copyFromAdaptiveVector(*x,*x_sg_blocks);
663  me_inargs.set_x_sg(x_sg_blocks);
664  }
665  if (x_dot != Teuchos::null) {
666  // copy to a poly orthog vector
667  adaptMngr->copyFromAdaptiveVector(*x_dot,*x_dot_sg_blocks);
668  me_inargs.set_x_dot_sg(x_dot_sg_blocks);
669  }
670  if (me_inargs.supports(IN_ARG_alpha))
671  me_inargs.set_alpha(inArgs.get_alpha());
672  if (me_inargs.supports(IN_ARG_beta))
673  me_inargs.set_beta(inArgs.get_beta());
674  if (me_inargs.supports(IN_ARG_t))
675  me_inargs.set_t(inArgs.get_t());
676  if (me_inargs.supports(IN_ARG_sg_basis)) {
677  if (inArgs.get_sg_basis() != Teuchos::null)
678  me_inargs.set_sg_basis(inArgs.get_sg_basis());
679  else
680  me_inargs.set_sg_basis(sg_basis);
681  }
682  if (me_inargs.supports(IN_ARG_sg_quadrature)) {
683  if (inArgs.get_sg_quadrature() != Teuchos::null)
684  me_inargs.set_sg_quadrature(inArgs.get_sg_quadrature());
685  else
686  me_inargs.set_sg_quadrature(sg_quad);
687  }
688  if (me_inargs.supports(IN_ARG_sg_expansion)) {
689  if (inArgs.get_sg_expansion() != Teuchos::null)
690  me_inargs.set_sg_expansion(inArgs.get_sg_expansion());
691  else
692  me_inargs.set_sg_expansion(sg_exp);
693  }
694 
695  // Pass parameters
696  for (int i=0; i<num_p; i++)
697  me_inargs.set_p(i, inArgs.get_p(i));
698  for (int i=0; i<num_p_sg; i++) {
699  Teuchos::RCP<const Epetra_Vector> p = inArgs.get_p(i+num_p);
700 
701  // We always need to pass in the SG parameters, so just use
702  // the initial parameters if it is NULL
703  if (p == Teuchos::null)
704  p = sg_p_init[i]->getBlockVector();
705 
706  // Convert block p to SG polynomial
708  create_p_sg(sg_p_index_map[i], View, p.get());
709  me_inargs.set_p_sg(sg_p_index_map[i], p_sg);
710  }
711 
712  // Create underlying outargs
713  OutArgs me_outargs = me->createOutArgs();
714 
715  // f
716  if (f_out != Teuchos::null) {
717  me_outargs.set_f_sg(f_sg_blocks);
718  if (eval_W_with_f)
719  me_outargs.set_W_sg(W_sg_blocks);
720  }
721 
722  // W
723  if (W_out != Teuchos::null && !eval_W_with_f )
724  me_outargs.set_W_sg(W_sg_blocks);
725 
726  // df/dp
727  for (int i=0; i<num_p; i++) {
728  if (!outArgs.supports(OUT_ARG_DfDp, i).none()) {
729  Derivative dfdp = outArgs.get_DfDp(i);
730  if (dfdp.getMultiVector() != Teuchos::null) {
732  if (dfdp.getMultiVectorOrientation() == DERIV_MV_BY_COL)
733  dfdp_sg =
735  sg_basis, overlapped_stoch_row_map,
736  me->get_f_map(), sg_comm,
737  me->get_p_map(i)->NumMyElements()));
738  else if (dfdp.getMultiVectorOrientation() == DERIV_TRANS_MV_BY_ROW)
739  dfdp_sg =
741  sg_basis, overlapped_stoch_row_map,
742  me->get_p_map(i), sg_comm,
743  me->get_f_map()->NumMyElements()));
744  me_outargs.set_DfDp_sg(i,
745  SGDerivative(dfdp_sg,
746  dfdp.getMultiVectorOrientation()));
747  }
748  TEUCHOS_TEST_FOR_EXCEPTION(dfdp.getLinearOp() != Teuchos::null, std::logic_error,
749  "Error! Stokhos::SGModelEvaluator_Adaptive::evalModel " <<
750  "cannot handle operator form of df/dp!");
751  }
752  }
753 
754  // Responses (g, dg/dx, dg/dp, ...)
755  for (int i=0; i<num_g_sg; i++) {
756  int ii = sg_g_index_map[i];
757 
758  // g
759  Teuchos::RCP<Epetra_Vector> g = outArgs.get_g(i);
760  if (g != Teuchos::null) {
762  create_g_sg(sg_g_index_map[i], View, g.get());
763  me_outargs.set_g_sg(ii, g_sg);
764  }
765 
766  // dg/dxdot
767  if (outArgs.supports(OUT_ARG_DgDx_dot, i).supports(DERIV_LINEAR_OP)) {
768  Derivative dgdx_dot = outArgs.get_DgDx_dot(i);
769  if (dgdx_dot.getLinearOp() != Teuchos::null) {
771  Teuchos::rcp_dynamic_cast<Stokhos::SGOperator>(
772  dgdx_dot.getLinearOp(), true);
774  op->getSGPolynomial();
775  if (me_outargs.supports(OUT_ARG_DgDx, ii).supports(DERIV_LINEAR_OP))
776  me_outargs.set_DgDx_dot_sg(ii, sg_blocks);
777  else {
778  for (unsigned int k=0; k<num_sg_blocks; k++) {
780  Teuchos::rcp_dynamic_cast<Stokhos::EpetraMultiVectorOperator>(
781  sg_blocks->getCoeffPtr(k), true)->getMultiVector();
782  dgdx_dot_sg_blocks[i]->setCoeffPtr(k, mv);
783  }
784  if (me_outargs.supports(OUT_ARG_DgDx_dot_sg, ii).supports(DERIV_MV_BY_COL))
785  me_outargs.set_DgDx_dot_sg(ii, SGDerivative(dgdx_dot_sg_blocks[i],
786  DERIV_MV_BY_COL));
787  else
788  me_outargs.set_DgDx_dot_sg(ii, SGDerivative(dgdx_dot_sg_blocks[i],
789  DERIV_TRANS_MV_BY_ROW));
790  }
791  }
792  TEUCHOS_TEST_FOR_EXCEPTION(dgdx_dot.getLinearOp() == Teuchos::null &&
793  dgdx_dot.isEmpty() == false,
794  std::logic_error,
795  "Error! Stokhos::SGModelEvaluator_Adaptive::evalModel: " <<
796  "Operator form of dg/dxdot is required!");
797  }
798 
799  // dg/dx
800  if (outArgs.supports(OUT_ARG_DgDx, i).supports(DERIV_LINEAR_OP)) {
801  Derivative dgdx = outArgs.get_DgDx(i);
802  if (dgdx.getLinearOp() != Teuchos::null) {
804  Teuchos::rcp_dynamic_cast<Stokhos::SGOperator>(
805  dgdx.getLinearOp(), true);
807  op->getSGPolynomial();
808  if (me_outargs.supports(OUT_ARG_DgDx, ii).supports(DERIV_LINEAR_OP))
809  me_outargs.set_DgDx_sg(ii, sg_blocks);
810  else {
811  for (unsigned int k=0; k<num_sg_blocks; k++) {
813  Teuchos::rcp_dynamic_cast<Stokhos::EpetraMultiVectorOperator>(
814  sg_blocks->getCoeffPtr(k), true)->getMultiVector();
815  dgdx_sg_blocks[i]->setCoeffPtr(k, mv);
816  }
817  if (me_outargs.supports(OUT_ARG_DgDx_sg, ii).supports(DERIV_MV_BY_COL))
818  me_outargs.set_DgDx_sg(ii, SGDerivative(dgdx_sg_blocks[i],
819  DERIV_MV_BY_COL));
820  else
821  me_outargs.set_DgDx_sg(ii, SGDerivative(dgdx_sg_blocks[i],
822  DERIV_TRANS_MV_BY_ROW));
823  }
824  }
825  TEUCHOS_TEST_FOR_EXCEPTION(dgdx.getLinearOp() == Teuchos::null &&
826  dgdx.isEmpty() == false,
827  std::logic_error,
828  "Error! Stokhos::SGModelEvaluator_Adaptive::evalModel: " <<
829  "Operator form of dg/dxdot is required!");
830  }
831 
832  // dg/dp
833  // Rembember, no derivatives w.r.t. sg parameters
834  for (int j=0; j<num_p; j++) {
835  if (!outArgs.supports(OUT_ARG_DgDp, i, j).none()) {
836  Derivative dgdp = outArgs.get_DgDp(i,j);
837  if (dgdp.getMultiVector() != Teuchos::null) {
839  if (dgdp.getMultiVectorOrientation() == DERIV_MV_BY_COL)
840  dgdp_sg =
842  sg_basis, overlapped_stoch_row_map,
843  me->get_g_map(ii), sg_g_map[i], sg_comm,
844  View, *(dgdp.getMultiVector())));
845  else if (dgdp.getMultiVectorOrientation() == DERIV_TRANS_MV_BY_ROW) {
847  Teuchos::rcp(&(dgdp.getMultiVector()->Map()),false);
848  dgdp_sg =
850  sg_basis, overlapped_stoch_row_map,
851  me->get_p_map(j), product_map, sg_comm,
852  View, *(dgdp.getMultiVector())));
853  }
854  me_outargs.set_DgDp_sg(ii, j,
855  SGDerivative(dgdp_sg,
856  dgdp.getMultiVectorOrientation()));
857  }
858  TEUCHOS_TEST_FOR_EXCEPTION(dgdp.getLinearOp() != Teuchos::null,
859  std::logic_error,
860  "Error! Stokhos::SGModelEvaluator_Adaptive::evalModel " <<
861  "cannot handle operator form of dg/dp!");
862  }
863  }
864 
865  }
866 
867  // Compute the functions
868  me->evalModel(me_inargs, me_outargs);
869 
870  // Copy block SG components for W
871  if ((W_out != Teuchos::null || (eval_W_with_f && f_out != Teuchos::null)) ) {
872 
874  if (W_out != Teuchos::null)
875  W = W_out;
876  else
877  W = my_W;
878 
879  Teuchos::RCP<Epetra_CrsMatrix> W_sg = Teuchos::rcp_dynamic_cast<Epetra_CrsMatrix>(W, true);
880  adaptMngr->setupOperator(*W_sg,*Cijk,*W_sg_blocks);
881  }
882 
883  // f
884  if (f_out!=Teuchos::null){
885  if (!scaleOP)
886  for (int i=0; i<sg_basis->size(); i++)
887  (*f_sg_blocks)[i].Scale(sg_basis->norm_squared(i));
888 
889  adaptMngr->copyToAdaptiveVector(*f_sg_blocks,*f_out);
890  }
891 
892  // df/dp
893  for (int i=0; i<num_p; i++) {
894  if (!outArgs.supports(OUT_ARG_DfDp, i).none()) {
895  Derivative dfdp = outArgs.get_DfDp(i);
896  SGDerivative dfdp_sg = me_outargs.get_DfDp_sg(i);
897  if (dfdp.getMultiVector() != Teuchos::null) {
898 
899  dfdp.getMultiVector()->Export(
900  *(dfdp_sg.getMultiVector()->getBlockMultiVector()),
901  *adapted_overlapped_f_exporter, Insert);
902  }
903  }
904  }
905 }
906 
907 void
909  const Stokhos::EpetraVectorOrthogPoly& x_sg_in)
910 {
911  *sg_x_init = x_sg_in;
912 }
913 
916 {
917  return sg_x_init;
918 }
919 
920 void
922  int i, const Stokhos::EpetraVectorOrthogPoly& p_sg_in)
923 {
924  *sg_p_init[i] = p_sg_in;
925 }
926 
929 {
930  return sg_p_init[l];
931 }
932 
935 {
936  return sg_p_index_map;
937 }
938 
941 {
942  return sg_g_index_map;
943 }
944 
947 {
949  for (int i=0; i<num_g; i++)
950  base_maps[i] = me->get_g_map(i);
951  return base_maps;
952  }
953 
956 {
957  return overlapped_stoch_row_map;
958 }
959 
962 {
963  return adapted_overlapped_x_map;
964 }
965 
968 {
969  return adapted_overlapped_x_importer;
970 }
971 
974  const Epetra_Vector* v) const
975 {
977  if (v == NULL)
979  sg_basis, stoch_row_map, x_map, sg_x_map, sg_comm));
980  else
982  sg_basis, stoch_row_map, x_map, sg_x_map, sg_comm,
983  CV, *v));
984  return sg_x;
985 }
986 
989  const Epetra_Vector* v) const
990 {
992  if (v == NULL)
994  sg_basis, overlapped_stoch_row_map, x_map,
995  sg_overlapped_x_map, sg_comm));
996  else
998  sg_basis, overlapped_stoch_row_map, x_map,
999  sg_overlapped_x_map, sg_comm, CV, *v));
1000  return sg_x;
1001 }
1002 
1005  const Epetra_MultiVector* v) const
1006 {
1008  if (v == NULL)
1010  sg_basis, stoch_row_map, x_map, sg_x_map, sg_comm,
1011  num_vecs));
1012  else
1014  sg_basis, stoch_row_map, x_map, sg_x_map, sg_comm,
1015  CV, *v));
1016  return sg_x;
1017 }
1018 
1021  int num_vecs,
1022  Epetra_DataAccess CV,
1023  const Epetra_MultiVector* v) const
1024 {
1026  if (v == NULL)
1028  sg_basis, overlapped_stoch_row_map, x_map,
1029  sg_overlapped_x_map, sg_comm, num_vecs));
1030  else
1032  sg_basis, overlapped_stoch_row_map, x_map,
1033  sg_overlapped_x_map, sg_comm, CV, *v));
1034  return sg_x;
1035 }
1036 
1039  const Epetra_Vector* v) const
1040 {
1042  Teuchos::Array<int>::const_iterator it = std::find(sg_p_index_map.begin(),
1043  sg_p_index_map.end(),
1044  l);
1045  TEUCHOS_TEST_FOR_EXCEPTION(it == sg_p_index_map.end(), std::logic_error,
1046  "Error! Invalid p map index " << l);
1047  int ll = it - sg_p_index_map.begin();
1048  if (v == NULL)
1050  sg_basis, overlapped_stoch_p_map, me->get_p_map(l),
1051  sg_p_map[ll], sg_comm));
1052  else
1054  sg_basis, overlapped_stoch_p_map, me->get_p_map(l),
1055  sg_p_map[ll], sg_comm, CV, *v));
1056  return sg_p;
1057 }
1058 
1061  const Epetra_Vector* v) const
1062 {
1064  if (v == NULL)
1066  sg_basis, stoch_row_map, f_map, sg_f_map, sg_comm));
1067  else
1069  sg_basis, stoch_row_map, f_map, sg_f_map, sg_comm,
1070  CV, *v));
1071  return sg_f;
1072 }
1073 
1076  const Epetra_Vector* v) const
1077 {
1079  if (v == NULL)
1081  sg_basis, overlapped_stoch_row_map, f_map,
1082  sg_overlapped_f_map, sg_comm));
1083  else
1085  sg_basis, overlapped_stoch_row_map, f_map,
1086  sg_overlapped_f_map, sg_comm, CV, *v));
1087  return sg_f;
1088 }
1089 
1092  int num_vecs,
1093  Epetra_DataAccess CV,
1094  const Epetra_MultiVector* v) const
1095 {
1097  if (v == NULL)
1099  sg_basis, stoch_row_map, f_map, sg_f_map, sg_comm,
1100  num_vecs));
1101  else
1103  sg_basis, stoch_row_map, f_map, sg_f_map, sg_comm,
1104  CV, *v));
1105  return sg_f;
1106 }
1107 
1110  int num_vecs,
1111  Epetra_DataAccess CV,
1112  const Epetra_MultiVector* v) const
1113 {
1115  if (v == NULL)
1117  sg_basis, overlapped_stoch_row_map, f_map,
1118  sg_overlapped_f_map, sg_comm, num_vecs));
1119  else
1121  sg_basis, overlapped_stoch_row_map, f_map,
1122  sg_overlapped_f_map, sg_comm, CV, *v));
1123  return sg_f;
1124 }
1125 
1128  const Epetra_Vector* v) const
1129 {
1131  Teuchos::Array<int>::const_iterator it = std::find(sg_g_index_map.begin(),
1132  sg_g_index_map.end(),
1133  l);
1134  TEUCHOS_TEST_FOR_EXCEPTION(it == sg_g_index_map.end(), std::logic_error,
1135  "Error! Invalid g map index " << l);
1136  int ll = it - sg_g_index_map.begin();
1137  if (v == NULL)
1139  sg_basis, overlapped_stoch_row_map,
1140  me->get_g_map(l),
1141  sg_g_map[ll], sg_comm));
1142  else
1144  sg_basis, overlapped_stoch_row_map,
1145  me->get_g_map(l),
1146  sg_g_map[ll], sg_comm, CV, *v));
1147  return sg_g;
1148 }
1149 
1152  Epetra_DataAccess CV,
1153  const Epetra_MultiVector* v) const
1154 {
1156  Teuchos::Array<int>::const_iterator it = std::find(sg_g_index_map.begin(),
1157  sg_g_index_map.end(),
1158  l);
1159  TEUCHOS_TEST_FOR_EXCEPTION(it == sg_g_index_map.end(), std::logic_error,
1160  "Error! Invalid g map index " << l);
1161  int ll = it - sg_g_index_map.begin();
1162  if (v == NULL)
1164  sg_basis, overlapped_stoch_row_map,
1165  me->get_g_map(l),
1166  sg_g_map[ll], sg_comm, num_vecs));
1167  else
1169  sg_basis, overlapped_stoch_row_map,
1170  me->get_g_map(l),
1171  sg_g_map[ll], sg_comm, CV, *v));
1172  return sg_g;
1173 }
Teuchos::RCP< const Epetra_Import > get_x_sg_importer() const
Return x sg importer.
Teuchos::Array< int > get_g_sg_map_indices() const
Get indices of SG responses.
Teuchos::RCP< const Epetra_Map > get_p_map(int l) const
Return parameter vector map.
bool eval_W_with_f
Whether to always evaluate W with f.
Teuchos::RCP< Stokhos::EpetraVectorOrthogPoly > x_sg_blocks
x stochastic Galerkin components
A container class storing an orthogonal polynomial whose coefficients are vectors, operators, or in general any type that would have an expensive copy constructor.
int num_g_sg
Number of stochastic response vectors.
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.
bool myGID(int i) const
Return whether global index i resides on this processor.
Teuchos::RCP< const Epetra_Map > sg_f_map
Block SG residual map.
Teuchos::RCP< const Teuchos::Array< std::string > > get_p_names(int l) const
Return array of parameter names.
Teuchos::RCP< const Epetra_BlockMap > stoch_row_map
Map for stochastic blocks.
Teuchos::RCP< const Stokhos::EpetraVectorOrthogPoly > get_p_sg_init(int l) const
Return initial SG parameters.
Teuchos::RCP< const Stokhos::EpetraSparse3Tensor > serialCijk
Serial Epetra Cijk for dgdx*.
T & get(ParameterList &l, const std::string &name)
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.
Teuchos::Array< Teuchos::RCP< Stokhos::EpetraMultiVectorOrthogPoly > > dgdx_sg_blocks
dg/dx stochastic Galerkin components
#define TEUCHOS_TEST_FOR_EXCEPTION(throw_exception_test, Exception, msg)
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.
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::RCP< Epetra_Export > adapted_overlapped_f_exporter
Exporter from SG-overlapped to SG maps.
Teuchos::Array< Teuchos::RCP< Teuchos::Array< std::string > > > sg_p_names
SG coefficient parameter names.
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< const Epetra_Map > sg_x_map
Block SG unknown map.
Teuchos::Array< Teuchos::RCP< Stokhos::EpetraMultiVectorOrthogPoly > > dgdx_dot_sg_blocks
dg/dxdot stochastic Galerkin components
void set_p_sg_init(int i, const Stokhos::EpetraVectorOrthogPoly &p_sg_in)
Set initial parameter polynomial.
T * get() const
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< const Stokhos::ParallelData > sg_parallel_data
Parallel SG data.
Teuchos::RCP< const Epetra_Map > adapted_overlapped_x_map
Adapated block SG overlapped unknown map.
Teuchos::RCP< const Epetra_BlockMap > get_x_sg_overlap_map() const
Return x sg overlap map.
virtual ordinal_type dimension() const =0
Return dimension of basis.
int num_g
Number of response vectors of underlying model evaluator.
Teuchos::RCP< Stokhos::EpetraVectorOrthogPoly > sg_x_init
SG initial x.
Teuchos::RCP< Epetra_Operator > create_W() const
Create W = alpha*M + beta*J matrix.
Teuchos::Array< Teuchos::RCP< const Epetra_Map > > sg_g_map
Block SG response map.
Teuchos::RCP< const Epetra_Map > get_g_map(int l) const
Return response map.
Teuchos::RCP< const Stokhos::EpetraSparse3Tensor > epetraCijk
Epetra Cijk.
Abstract base class for orthogonal polynomial-based expansions.
unsigned int num_W_blocks
Number of W stochastic blocks (may be smaller than num_sg_blocks)
Teuchos::Array< int > get_p_sg_map_indices() const
Get indices of SG parameters.
std::vector< Teuchos::RCP< const Stokhos::ProductBasis< int, double > > > sg_row_dof_basis
Teuchos::RCP< EpetraExt::ModelEvaluator > me
Underlying model evaluator.
An adaptor that supplies the operator interface to a multi-vector.
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::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.
bool supports_x
Whether we support x (and thus f and W)
int num_p
Number of parameter vectors of underlying model evaluator.
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::RCP< Stokhos::EpetraVectorOrthogPoly > create_p_sg(int l, Epetra_DataAccess CV=Copy, const Epetra_Vector *v=0) const
Create vector orthog poly using p map.
Teuchos::RCP< const Epetra_Map > sg_overlapped_x_map
Block SG overlapped unknown map.
virtual Teuchos::RCP< Stokhos::EpetraOperatorOrthogPoly > getSGPolynomial()=0
Get SG polynomial.
void evalModel(const InArgs &inArgs, const OutArgs &outArgs) const
Evaluate model on InArgs.
Teuchos::RCP< Stokhos::EpetraVectorOrthogPoly > x_dot_sg_blocks
x_dot stochastic Galerkin components
Teuchos::RCP< const Epetra_Map > adapted_overlapped_f_map
Adapted block SG overlapped residual map.
Teuchos::RCP< const Epetra_BlockMap > overlapped_stoch_row_map
Overlapped map for stochastic blocks (local map)
Teuchos::RCP< const Stokhos::EpetraVectorOrthogPoly > get_x_sg_init() const
Return initial SG x.
Teuchos::RCP< const Epetra_Map > adapted_f_map
Adapted block SG residual map.
iterator end()
std::vector< T >::const_iterator const_iterator
Teuchos::RCP< const Epetra_Map > adapted_x_map
Adapted lock SG unknown map.
Teuchos::Array< Teuchos::RCP< const Epetra_Map > > sg_p_map
Block SG parameter map.
Teuchos::Array< Teuchos::RCP< const Epetra_Map > > get_g_sg_base_maps() const
Get base maps of SG responses.
void push_back(const value_type &x)
Teuchos::RCP< Stokhos::EpetraOperatorOrthogPoly > W_sg_blocks
W stochastic Galerkin components.
Teuchos::RCP< Epetra_Import > adapted_overlapped_x_importer
Importer from SG to SG-overlapped maps.
Teuchos::RCP< const Epetra_BlockMap > get_overlap_stochastic_map() const
Return overlap stochastic map.
Teuchos::RCP< const Epetra_Map > get_x_map() const
Return solution vector map.
Teuchos::Array< int > sg_p_index_map
Index map between block-p and p_sg maps.
size_type size() const
int num_p_sg
Number of stochastic parameter vectors.
void set_x_sg_init(const Stokhos::EpetraVectorOrthogPoly &x_sg_in)
Set initial solution polynomial.
unsigned int num_sg_blocks
Number of stochastic blocks.
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< 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.
Teuchos::RCP< const Epetra_Vector > get_p_init(int l) const
Return initial parameters.
const Epetra_CrsGraph & Graph() const
Teuchos::Array< int > sg_g_index_map
Index map between block-g and g_sg maps.
int Scale(double ScalarConstant)
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.
Epetra_DataAccess
Teuchos::RCP< const EpetraExt::MultiComm > sg_comm
Parallel SG communicator.
Teuchos::RCP< Epetra_Vector > my_x
x pointer for evaluating preconditioner
Teuchos::RCP< const Stokhos::OrthogPolyBasis< int, double > > sg_basis
Stochastic Galerkin basis.
Teuchos::RCP< Teuchos::ParameterList > params
Algorithmic parameters.
unsigned int num_p_blocks
Number of p stochastic blocks (may be smaller than num_sg_blocks)
Teuchos::RCP< Stokhos::AdaptivityManager > adaptMngr
iterator begin()
Teuchos::RCP< const Epetra_Map > sg_overlapped_f_map
Block SG overlapped residual map.
ScalarType g(const Teuchos::Array< ScalarType > &x, const ScalarType &y)
Teuchos::RCP< Stokhos::EpetraVectorOrthogPoly > f_sg_blocks
f stochastic Galerkin components
Teuchos::RCP< const Epetra_Vector > get_x_init() const
Return initial solution.
Teuchos::RCP< const Epetra_BlockMap > overlapped_stoch_p_map
Overlapped map for p stochastic blocks (local map)
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.
SGModelEvaluator_Adaptive(const Teuchos::RCP< EpetraExt::ModelEvaluator > &me_, const Teuchos::RCP< Stokhos::AdaptivityManager > &am, 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_, bool onlyUseLinear_, int kExpOrder_, const Teuchos::RCP< Teuchos::ParameterList > &params_)
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< Teuchos::RCP< Stokhos::EpetraVectorOrthogPoly > > sg_p_init
SG initial p.
Teuchos::RCP< const Epetra_Map > f_map
Underlying residual map.