Stokhos Package Browser (Single Doxygen Collection)  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Stokhos_SGQuadModelEvaluator.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 
12 #include "Stokhos_Quadrature.hpp"
16 #include "Epetra_Map.h"
17 #include "Epetra_Vector.h"
18 #include "Teuchos_TimeMonitor.hpp"
19 #include "Teuchos_Assert.hpp"
20 
24  me(me_),
25  num_p(0),
26  num_g(0),
27  x_dot_qp(),
28  x_qp(),
29  p_qp(),
30  f_qp(),
31  W_qp(),
32  dfdp_qp(),
33  g_qp(),
34  dgdx_qp(),
35  dgdx_dot_qp(),
36  dgdp_qp()
37 {
38  // Create storage for x_dot, x, and p at a quad point
39  InArgs me_inargs = me->createInArgs();
40  num_p = me_inargs.Np();
41  if (me_inargs.supports(IN_ARG_x_dot))
42  x_dot_qp = Teuchos::rcp(new Epetra_Vector(*(me->get_x_map())));
43  if (me_inargs.supports(IN_ARG_x))
44  x_qp = Teuchos::rcp(new Epetra_Vector((*me->get_x_map())));
45  p_qp.resize(num_p);
46  for (int i=0; i<num_p; i++)
47  p_qp[i] = Teuchos::rcp(new Epetra_Vector(*(me->get_p_map(i))));
48 
49  // Create storage for f and W at a quad point
50  OutArgs me_outargs = me->createOutArgs();
51  num_g = me_outargs.Ng();
52 
53  // f
54  if (me_outargs.supports(OUT_ARG_f))
55  f_qp = Teuchos::rcp(new Epetra_Vector(*(me->get_f_map())));
56 
57  // W
58  if (me_outargs.supports(OUT_ARG_W))
59  W_qp = me->create_W();
60 
61  // df/dp
62  dfdp_qp.resize(num_p);
63  for (int i=0; i<num_p; i++)
64  if (me_outargs.supports(OUT_ARG_DfDp,i).supports(DERIV_MV_BY_COL))
65  dfdp_qp[i] = EpetraExt::ModelEvaluator::Derivative(
67  *(me->get_f_map()),
68  me->get_p_map(i)->NumGlobalElements())));
69  else if (me_outargs.supports(OUT_ARG_DfDp,i).supports(DERIV_TRANS_MV_BY_ROW))
70  dfdp_qp[i] = EpetraExt::ModelEvaluator::Derivative(
72  *(me->get_p_map(i)),
73  me->get_f_map()->NumGlobalElements())));
74  else if (me_outargs.supports(OUT_ARG_DfDp,i).supports(DERIV_LINEAR_OP))
75  dfdp_qp[i] = EpetraExt::ModelEvaluator::Derivative(
76  me->create_DfDp_op(i));
77 
78  g_qp.resize(num_g);
81  dgdp_qp.resize(num_g);
82  for (int i=0; i<num_g; i++) {
83 
84  // g
85  g_qp[i] =
86  Teuchos::rcp(new Epetra_Vector(*(me->get_g_map(i))));
87 
88  // dg/dx
89  if (me_outargs.supports(OUT_ARG_DgDx, i).supports(DERIV_TRANS_MV_BY_ROW))
90  dgdx_qp[i] = EpetraExt::ModelEvaluator::Derivative(
92  *(me->get_x_map()),
93  me->get_g_map(i)->NumGlobalElements())));
94  else if (me_outargs.supports(OUT_ARG_DgDx, i).supports(DERIV_MV_BY_COL))
95  dgdx_qp[i] = EpetraExt::ModelEvaluator::Derivative(
97  *(me->get_g_map(i)),
98  me->get_x_map()->NumGlobalElements())));
99  else if (me_outargs.supports(OUT_ARG_DgDx, i).supports(DERIV_LINEAR_OP))
100  dgdx_qp[i] = EpetraExt::ModelEvaluator::Derivative(
101  me->create_DgDx_op(i));
102 
103  // dg/dx_dot
104  if (me_outargs.supports(OUT_ARG_DgDx_dot, i).supports(DERIV_TRANS_MV_BY_ROW))
105  dgdx_dot_qp[i] = EpetraExt::ModelEvaluator::Derivative(
107  *(me->get_x_map()),
108  me->get_g_map(i)->NumGlobalElements())));
109  else if (me_outargs.supports(OUT_ARG_DgDx_dot, i).supports(DERIV_MV_BY_COL))
110  dgdx_dot_qp[i] = EpetraExt::ModelEvaluator::Derivative(
112  *(me->get_g_map(i)),
113  me->get_x_map()->NumGlobalElements())));
114  else if (me_outargs.supports(OUT_ARG_DgDx_dot, i).supports(DERIV_LINEAR_OP))
115  dgdx_dot_qp[i] = EpetraExt::ModelEvaluator::Derivative(
116  me->create_DgDx_dot_op(i));
117 
118  // dg/dp
119  dgdp_qp[i].resize(num_p);
120  for (int j=0; j<num_p; j++)
121  if (me_outargs.supports(OUT_ARG_DgDp, i, j).supports(DERIV_TRANS_MV_BY_ROW))
122  dgdp_qp[i][j] = EpetraExt::ModelEvaluator::Derivative(
124  *(me->get_p_map(j)),
125  me->get_g_map(i)->NumGlobalElements())));
126  else if (me_outargs.supports(OUT_ARG_DgDp, i, j).supports(DERIV_MV_BY_COL))
127  dgdp_qp[i][j] = EpetraExt::ModelEvaluator::Derivative(
129  *(me->get_g_map(i)),
130  me->get_p_map(j)->NumGlobalElements())));
131  else if (me_outargs.supports(OUT_ARG_DgDp, i, j).supports(DERIV_LINEAR_OP))
132  dgdp_qp[i][j] = EpetraExt::ModelEvaluator::Derivative(
133  me->create_DgDp_op(i,j));
134  }
135 }
136 
137 // Overridden from EpetraExt::ModelEvaluator
138 
141 get_x_map() const
142 {
143  return me->get_x_map();
144 }
145 
148 get_f_map() const
149 {
150  return me->get_f_map();
151 }
152 
155 get_p_map(int l) const
156 {
157  return me->get_p_map(l);
158 }
159 
162 get_g_map(int l) const
163 {
164  return me->get_g_map(l);
165 }
166 
169 get_p_names(int l) const
170 {
171  return me->get_p_names(l);
172 }
173 
176 get_x_init() const
177 {
178  return me->get_x_init();
179 }
180 
183 get_p_init(int l) const
184 {
185  return me->get_p_init(l);
186 }
187 
190 create_W() const
191 {
192  return me->create_W();
193 }
194 
195 EpetraExt::ModelEvaluator::InArgs
198 {
199  InArgsSetup inArgs;
200  InArgs me_inargs = me->createInArgs();
201 
202  inArgs.setModelEvalDescription(this->description());
203  inArgs.set_Np(num_p);
204  inArgs.setSupports(IN_ARG_x_dot, me_inargs.supports(IN_ARG_x_dot));
205  inArgs.setSupports(IN_ARG_x, me_inargs.supports(IN_ARG_x));
206  inArgs.setSupports(IN_ARG_t, me_inargs.supports(IN_ARG_t));
207  inArgs.setSupports(IN_ARG_alpha, me_inargs.supports(IN_ARG_alpha));
208  inArgs.setSupports(IN_ARG_beta, me_inargs.supports(IN_ARG_beta));
209 
210  for (int i=0; i<num_p; i++)
211  inArgs.setSupports(IN_ARG_p_sg, i, true);
212  inArgs.setSupports(IN_ARG_x_sg, me_inargs.supports(IN_ARG_x));
213  inArgs.setSupports(IN_ARG_x_dot_sg, me_inargs.supports(IN_ARG_x_dot));
214  inArgs.setSupports(IN_ARG_sg_basis, true);
215  inArgs.setSupports(IN_ARG_sg_quadrature, true);
216 
217  return inArgs;
218 }
219 
220 EpetraExt::ModelEvaluator::OutArgs
223 {
224  OutArgsSetup outArgs;
225  OutArgs me_outargs = me->createOutArgs();
226 
227  outArgs.setModelEvalDescription(this->description());
228  outArgs.set_Np_Ng(num_p, num_g);
229  outArgs.setSupports(OUT_ARG_f, me_outargs.supports(OUT_ARG_f));
230  outArgs.setSupports(OUT_ARG_W, me_outargs.supports(OUT_ARG_W));
231  for (int j=0; j<num_p; j++)
232  outArgs.setSupports(OUT_ARG_DfDp, j,
233  me_outargs.supports(OUT_ARG_DfDp, j));
234  for (int i=0; i<num_g; i++) {
235  outArgs.setSupports(OUT_ARG_DgDx, i,
236  me_outargs.supports(OUT_ARG_DgDx, i));
237  outArgs.setSupports(OUT_ARG_DgDx_dot, i,
238  me_outargs.supports(OUT_ARG_DgDx_dot, i));
239  for (int j=0; j<num_p; j++)
240  outArgs.setSupports(OUT_ARG_DgDp, i, j,
241  me_outargs.supports(OUT_ARG_DgDp, i, j));
242  }
243 
244  outArgs.setSupports(OUT_ARG_f_sg, me_outargs.supports(OUT_ARG_f));
245  if (me_outargs.supports(OUT_ARG_W)) {
246  outArgs.set_W_properties(me_outargs.get_W_properties());
247  outArgs.setSupports(OUT_ARG_W_sg, true);
248  }
249  for (int j=0; j<num_p; j++) {
250  outArgs.setSupports(OUT_ARG_DfDp_sg, j,
251  me_outargs.supports(OUT_ARG_DfDp, j));
252  }
253  for (int i=0; i<num_g; i++) {
254  outArgs.setSupports(OUT_ARG_g_sg, i, true);
255  outArgs.setSupports(OUT_ARG_DgDx_sg, i,
256  me_outargs.supports(OUT_ARG_DgDx, i));
257  outArgs.setSupports(OUT_ARG_DgDx_dot_sg, i,
258  me_outargs.supports(OUT_ARG_DgDx_dot, i));
259  for (int j=0; j<num_p; j++)
260  outArgs.setSupports(OUT_ARG_DgDp_sg, i, j,
261  me_outargs.supports(OUT_ARG_DgDp, i, j));
262  }
263 
264  return outArgs;
265 }
266 
267 void
269 evalModel(const InArgs& inArgs, const OutArgs& outArgs) const
270 {
271  // Create underlying inargs
272  InArgs me_inargs = me->createInArgs();
273  if (me_inargs.supports(IN_ARG_x))
274  me_inargs.set_x(inArgs.get_x());
275  if (me_inargs.supports(IN_ARG_x_dot))
276  me_inargs.set_x_dot(inArgs.get_x_dot());
277  if (me_inargs.supports(IN_ARG_alpha))
278  me_inargs.set_alpha(inArgs.get_alpha());
279  if (me_inargs.supports(IN_ARG_beta))
280  me_inargs.set_beta(inArgs.get_beta());
281  if (me_inargs.supports(IN_ARG_t))
282  me_inargs.set_t(inArgs.get_t());
283  for (int i=0; i<num_p; i++)
284  me_inargs.set_p(i, inArgs.get_p(i));
285 
286  // Create underlying outargs
287  OutArgs me_outargs = me->createOutArgs();
288  if (me_outargs.supports(OUT_ARG_f))
289  me_outargs.set_f(outArgs.get_f());
290  if (me_outargs.supports(OUT_ARG_W))
291  me_outargs.set_W(outArgs.get_W());
292  for (int j=0; j<num_p; j++)
293  if (!outArgs.supports(OUT_ARG_DfDp, j).none())
294  me_outargs.set_DfDp(j, outArgs.get_DfDp(j));
295  for (int i=0; i<num_g; i++) {
296  me_outargs.set_g(i, outArgs.get_g(i));
297  if (!outArgs.supports(OUT_ARG_DgDx, i).none())
298  me_outargs.set_DgDx(i, outArgs.get_DgDx(i));
299  if (!outArgs.supports(OUT_ARG_DgDx_dot, i).none())
300  me_outargs.set_DgDx(i, outArgs.get_DgDx_dot(i));
301  for (int j=0; j<num_p; j++)
302  if (!outArgs.supports(OUT_ARG_DgDp, i, j).none())
303  me_outargs.set_DgDp(i, j, outArgs.get_DgDp(i,j));
304  }
305 
306  bool do_quad = false;
307  InArgs::sg_const_vector_t x_sg;
308  InArgs::sg_const_vector_t x_dot_sg;
310  OutArgs::sg_vector_t f_sg;
311  OutArgs::sg_operator_t W_sg;
312  Teuchos::Array<SGDerivative> dfdp_sg(num_p);
314  Teuchos::Array<SGDerivative> dgdx_sg(num_g);
315  Teuchos::Array<SGDerivative> dgdx_dot_sg(num_g);
317  TEUCHOS_TEST_FOR_EXCEPTION(inArgs.get_sg_basis() == Teuchos::null,
318  std::logic_error,
319  "Error! Stokhos::SGQuadModelEvaluator::evalModel(): " <<
320  "SG basis inArg cannot be null!");
321  TEUCHOS_TEST_FOR_EXCEPTION(inArgs.get_sg_quadrature() == Teuchos::null,
322  std::logic_error,
323  "Error! Stokhos::SGQuadModelEvaluator::evalModel(): " <<
324  "SG quadrature inArg cannot be null!");
326  inArgs.get_sg_basis();
328  inArgs.get_sg_quadrature();
329  if (inArgs.supports(IN_ARG_x_sg)) {
330  x_sg = inArgs.get_x_sg();
331  if (x_sg != Teuchos::null) {
332  do_quad = true;
333  }
334  }
335  if (inArgs.supports(IN_ARG_x_dot_sg)) {
336  x_dot_sg = inArgs.get_x_dot_sg();
337  if (x_dot_sg != Teuchos::null) {
338  do_quad = true;
339  }
340  }
341  for (int i=0; i<num_p; i++) {
342  p_sg[i] = inArgs.get_p_sg(i);
343  if (p_sg[i] != Teuchos::null) {
344  do_quad = true;
345  }
346  }
347  if (outArgs.supports(OUT_ARG_f_sg)) {
348  f_sg = outArgs.get_f_sg();
349  if (f_sg != Teuchos::null)
350  f_sg->init(0.0);
351  }
352  if (outArgs.supports(OUT_ARG_W_sg)) {
353  W_sg = outArgs.get_W_sg();
354  if (W_sg != Teuchos::null)
355  W_sg->init(0.0);
356  }
357  for (int i=0; i<num_p; i++) {
358  if (!outArgs.supports(OUT_ARG_DfDp_sg, i).none()) {
359  dfdp_sg[i] = outArgs.get_DfDp_sg(i);
360  if (dfdp_sg[i].getMultiVector() != Teuchos::null)
361  dfdp_sg[i].getMultiVector()->init(0.0);
362  else if (dfdp_sg[i].getLinearOp() != Teuchos::null)
363  dfdp_sg[i].getLinearOp()->init(0.0);
364  }
365  }
366 
367  for (int i=0; i<num_g; i++) {
368  g_sg[i] = outArgs.get_g_sg(i);
369  if (g_sg[i] != Teuchos::null)
370  g_sg[i]->init(0.0);
371 
372  if (!outArgs.supports(OUT_ARG_DgDx_sg, i).none()) {
373  dgdx_sg[i] = outArgs.get_DgDx_sg(i);
374  if (dgdx_sg[i].getMultiVector() != Teuchos::null)
375  dgdx_sg[i].getMultiVector()->init(0.0);
376  else if (dgdx_sg[i].getLinearOp() != Teuchos::null)
377  dgdx_sg[i].getLinearOp()->init(0.0);
378  }
379 
380  if (!outArgs.supports(OUT_ARG_DgDx_dot_sg, i).none()) {
381  dgdx_dot_sg[i] = outArgs.get_DgDx_dot_sg(i);
382  if (dgdx_dot_sg[i].getMultiVector() != Teuchos::null)
383  dgdx_dot_sg[i].getMultiVector()->init(0.0);
384  else if (dgdx_dot_sg[i].getLinearOp() != Teuchos::null)
385  dgdx_dot_sg[i].getLinearOp()->init(0.0);
386  }
387 
388  dgdp_sg[i].resize(num_p);
389  for (int j=0; j<num_p; j++) {
390  if (!outArgs.supports(OUT_ARG_DgDp_sg, i, j).none()) {
391  dgdp_sg[i][j] = outArgs.get_DgDp_sg(i,j);
392  if (dgdp_sg[i][j].getMultiVector() != Teuchos::null)
393  dgdp_sg[i][j].getMultiVector()->init(0.0);
394  else if (dgdp_sg[i][j].getLinearOp() != Teuchos::null)
395  dgdp_sg[i][j].getLinearOp()->init(0.0);
396  }
397  }
398  }
399 
400  if (do_quad) {
401  // Get quadrature data
402  const Teuchos::Array< Teuchos::Array<double> >& quad_points =
403  quad->getQuadPoints();
404  const Teuchos::Array<double>& quad_weights =
405  quad->getQuadWeights();
406  const Teuchos::Array< Teuchos::Array<double> > & quad_values =
407  quad->getBasisAtQuadPoints();
408  const Teuchos::Array<double>& basis_norms = basis->norm_squared();
409 
410  // Perform integrations
411  for (int qp=0; qp<quad_points.size(); qp++) {
412 
413  // StieltjesPCEBasis can introduce quadrature points with zero weight
414  // Don't do those evaluations, since the model might not like the
415  // quadrature points (i.e., zero)
416  if (quad_weights[qp] == 0.0)
417  continue;
418 
419  {
420 #ifdef STOKHOS_TEUCHOS_TIME_MONITOR
421  TEUCHOS_FUNC_TIME_MONITOR_DIFF("Stokhos: SGQuadModelEvaluator -- Polynomial Evaluation",
422  PolyEvaluation);
423 #endif
424 
425  // Evaluate inputs at quadrature points
426  if (x_sg != Teuchos::null) {
427 #ifdef STOKHOS_TEUCHOS_TIME_MONITOR
428  TEUCHOS_FUNC_TIME_MONITOR("Stokhos: SGQuadModelEvaluator -- X Evaluation");
429 #endif
430  x_sg->evaluate(quad_values[qp], *x_qp);
431  me_inargs.set_x(x_qp);
432  }
433  if (x_dot_sg != Teuchos::null) {
434 #ifdef STOKHOS_TEUCHOS_TIME_MONITOR
435  TEUCHOS_FUNC_TIME_MONITOR("Stokhos: SGQuadModelEvaluator -- X_dot Evaluation");
436 #endif
437  x_dot_sg->evaluate(quad_values[qp], *x_dot_qp);
438  me_inargs.set_x_dot(x_qp);
439  }
440  for (int i=0; i<num_p; i++) {
441  if (p_sg[i] != Teuchos::null) {
442 #ifdef STOKHOS_TEUCHOS_TIME_MONITOR
443  TEUCHOS_FUNC_TIME_MONITOR("Stokhos: SGQuadModelEvaluator -- P Evaluation");
444 #endif
445  p_sg[i]->evaluate(quad_values[qp], *(p_qp[i]));
446  me_inargs.set_p(i, p_qp[i]);
447  }
448  }
449  if (f_sg != Teuchos::null)
450  me_outargs.set_f(f_qp);
451  if (W_sg != Teuchos::null)
452  me_outargs.set_W(W_qp);
453  for (int i=0; i<num_p; i++) {
454  if (!dfdp_sg[i].isEmpty())
455  me_outargs.set_DfDp(i, dfdp_qp[i]);
456  }
457  for (int i=0; i<num_g; i++) {
458  if (g_sg[i] != Teuchos::null)
459  me_outargs.set_g(i, g_qp[i]);
460  if (!dgdx_dot_sg[i].isEmpty())
461  me_outargs.set_DgDx_dot(i, dgdx_dot_qp[i]);
462  if (!dgdx_sg[i].isEmpty())
463  me_outargs.set_DgDx(i, dgdx_qp[i]);
464  for (int j=0; j<num_p; j++)
465  if (!dgdp_sg[i][j].isEmpty())
466  me_outargs.set_DgDp(i, j, dgdp_qp[i][j]);
467  }
468 
469  }
470 
471  {
472 #ifdef STOKHOS_TEUCHOS_TIME_MONITOR
473  TEUCHOS_FUNC_TIME_MONITOR("Stokhos: SGQuadModelEvaluator -- Model Evaluation");
474 #endif
475 
476  // Evaluate model at quadrature points
477  me->evalModel(me_inargs, me_outargs);
478 
479  }
480 
481  {
482 #ifdef STOKHOS_TEUCHOS_TIME_MONITOR
484  "Stokhos: SGQuadModelEvaluator -- Polynomial Integration", Integration);
485 #endif
486 
487  // Sum in results
488  if (f_sg != Teuchos::null) {
489 #ifdef STOKHOS_TEUCHOS_TIME_MONITOR
490  TEUCHOS_FUNC_TIME_MONITOR("Stokhos: SGQuadModelEvaluator -- F Integration");
491 #endif
492  f_sg->sumIntoAllTerms(quad_weights[qp], quad_values[qp], basis_norms,
493  *f_qp);
494  }
495  if (W_sg != Teuchos::null) {
496 #ifdef STOKHOS_TEUCHOS_TIME_MONITOR
497  TEUCHOS_FUNC_TIME_MONITOR("Stokhos: SGQuadModelEvaluator -- W Integration");
498 #endif
499  W_sg->sumIntoAllTerms(quad_weights[qp], quad_values[qp], basis_norms,
500  *W_qp);
501  }
502  for (int j=0; j<num_p; j++) {
503  if (!dfdp_sg[j].isEmpty()) {
504 #ifdef STOKHOS_TEUCHOS_TIME_MONITOR
506  "Stokhos: SGQuadModelEvaluator -- df/dp Integration");
507 #endif
508  if (dfdp_sg[j].getMultiVector() != Teuchos::null) {
509  dfdp_sg[j].getMultiVector()->sumIntoAllTerms(
510  quad_weights[qp], quad_values[qp], basis_norms,
511  *(dfdp_qp[j].getMultiVector()));
512  }
513  else if (dfdp_sg[j].getLinearOp() != Teuchos::null) {
514  dfdp_sg[j].getLinearOp()->sumIntoAllTerms(
515  quad_weights[qp], quad_values[qp], basis_norms,
516  *(dfdp_qp[j].getLinearOp()));
517  }
518  }
519  }
520  for (int i=0; i<num_g; i++) {
521  if (g_sg[i] != Teuchos::null) {
522 #ifdef STOKHOS_TEUCHOS_TIME_MONITOR
523  TEUCHOS_FUNC_TIME_MONITOR("Stokhos: SGQuadModelEvaluator -- G Integration");
524 #endif
525  g_sg[i]->sumIntoAllTerms(quad_weights[qp], quad_values[qp],
526  basis_norms, *g_qp[i]);
527  }
528  if (!dgdx_dot_sg[i].isEmpty()) {
529 #ifdef STOKHOS_TEUCHOS_TIME_MONITOR
531  "Stokhos: SGQuadModelEvaluator -- dg/dx_dot Integration");
532 #endif
533  if (dgdx_dot_sg[i].getMultiVector() != Teuchos::null) {
534  dgdx_dot_sg[i].getMultiVector()->sumIntoAllTerms(
535  quad_weights[qp], quad_values[qp], basis_norms,
536  *(dgdx_dot_qp[i].getMultiVector()));
537  }
538  else if (dgdx_dot_sg[i].getLinearOp() != Teuchos::null) {
539  dgdx_dot_sg[i].getLinearOp()->sumIntoAllTerms(
540  quad_weights[qp], quad_values[qp], basis_norms,
541  *(dgdx_dot_qp[i].getLinearOp()));
542  }
543  }
544  if (!dgdx_sg[i].isEmpty()) {
545 #ifdef STOKHOS_TEUCHOS_TIME_MONITOR
547  "Stokhos: SGQuadModelEvaluator -- dg/dx Integration");
548 #endif
549  if (dgdx_sg[i].getMultiVector() != Teuchos::null) {
550  dgdx_sg[i].getMultiVector()->sumIntoAllTerms(
551  quad_weights[qp], quad_values[qp], basis_norms,
552  *(dgdx_qp[i].getMultiVector()));
553  }
554  else if (dgdx_sg[i].getLinearOp() != Teuchos::null) {
555  dgdx_sg[i].getLinearOp()->sumIntoAllTerms(
556  quad_weights[qp], quad_values[qp], basis_norms,
557  *(dgdx_qp[i].getLinearOp()));
558  }
559  }
560  for (int j=0; j<num_p; j++) {
561  if (!dgdp_sg[i][j].isEmpty()) {
562 #ifdef STOKHOS_TEUCHOS_TIME_MONITOR
564  "Stokhos: SGQuadModelEvaluator -- dg/dp Integration");
565 #endif
566  if (dgdp_sg[i][j].getMultiVector() != Teuchos::null) {
567  dgdp_sg[i][j].getMultiVector()->sumIntoAllTerms(
568  quad_weights[qp], quad_values[qp], basis_norms,
569  *(dgdp_qp[i][j].getMultiVector()));
570  }
571  else if (dgdp_sg[i][j].getLinearOp() != Teuchos::null) {
572  dgdp_sg[i][j].getLinearOp()->sumIntoAllTerms(
573  quad_weights[qp], quad_values[qp], basis_norms,
574  *(dgdp_qp[i][j].getLinearOp()));
575  }
576  }
577  }
578  }
579 
580  }
581  }
582  }
583  else {
584  // Compute the non-SG functions
585  me->evalModel(me_inargs, me_outargs);
586  }
587 }
int num_p
Number of parameter vectors.
#define TEUCHOS_FUNC_TIME_MONITOR(FUNCNAME)
Teuchos::RCP< const Epetra_Map > get_x_map() const
Return solution vector map.
int num_g
Number of response vectors.
Teuchos::RCP< const Epetra_Vector > get_x_init() const
Return initial solution.
#define TEUCHOS_TEST_FOR_EXCEPTION(throw_exception_test, Exception, msg)
Teuchos::Array< EpetraExt::ModelEvaluator::Derivative > dgdx_dot_qp
Response derivative.
Teuchos::Array< EpetraExt::ModelEvaluator::Derivative > dgdx_qp
Response derivative.
Teuchos::RCP< const Teuchos::Array< std::string > > get_p_names(int l) const
Return array of parameter names.
InArgs createInArgs() const
Create InArgs.
Teuchos::Array< Teuchos::RCP< Epetra_Vector > > g_qp
Response vectors.
TEUCHOS_DEPRECATED RCP< T > rcp(T *p, Dealloc_T dealloc, bool owns_mem)
Teuchos::Array< Teuchos::RCP< Epetra_Vector > > p_qp
Parameter vectors.
Teuchos::RCP< const Epetra_Map > get_p_map(int l) const
Return parameter vector map.
Teuchos::Array< EpetraExt::ModelEvaluator::Derivative > dfdp_qp
Residual derivatives.
Teuchos::RCP< Epetra_Vector > x_dot_qp
Time derivative vector.
void resize(size_type new_size, const value_type &x=value_type())
Teuchos::RCP< Epetra_Operator > create_W() const
Create W = alpha*M + beta*J matrix.
Teuchos::RCP< Epetra_Vector > f_qp
Residual vector.
Teuchos::RCP< Epetra_Vector > x_qp
Solution vector.
Teuchos::RCP< Epetra_Operator > W_qp
W operator.
size_type size() const
Teuchos::Array< Teuchos::Array< EpetraExt::ModelEvaluator::Derivative > > dgdp_qp
Response sensitivities.
Teuchos::RCP< EpetraExt::ModelEvaluator > me
Underlying model evaluator.
Teuchos::RCP< const Epetra_Vector > get_p_init(int l) const
Return initial parameters.
Teuchos::RCP< const Epetra_Map > get_f_map() const
Return residual vector map.
OutArgs createOutArgs() const
Create OutArgs.
SGQuadModelEvaluator(const Teuchos::RCP< EpetraExt::ModelEvaluator > &me)
void evalModel(const InArgs &inArgs, const OutArgs &outArgs) const
Evaluate model on InArgs.
#define TEUCHOS_FUNC_TIME_MONITOR_DIFF(FUNCNAME, DIFF)
Teuchos::RCP< const Epetra_Map > get_g_map(int l) const
Return observation vector map.