Thyra  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Thyra_ModelEvaluatorHelpers.hpp
1 // @HEADER
2 // *****************************************************************************
3 // Thyra: Interfaces and Support for Abstract Numerical Algorithms
4 //
5 // Copyright 2004 NTESS and the Thyra contributors.
6 // SPDX-License-Identifier: BSD-3-Clause
7 // *****************************************************************************
8 // @HEADER
9 
10 #ifndef THYRA_MODEL_EVALUATOR_HELPERS_HPP
11 #define THYRA_MODEL_EVALUATOR_HELPERS_HPP
12 
13 
14 #include "Thyra_ModelEvaluator.hpp"
15 
16 
17 namespace Thyra {
18 
19 
27 template<class Scalar>
28 RCP<ModelEvaluatorBase::InArgs<Scalar> >
29 clone( const ModelEvaluatorBase::InArgs<Scalar> &inArgs );
30 
31 
33 template<class Scalar>
34 ModelEvaluatorBase::Derivative<Scalar>
35 derivativeGradient(
36  const RCP<MultiVectorBase<Scalar> > &grad
37  );
38 
39 
41 template<class Scalar>
42 ModelEvaluatorBase::DerivativeMultiVector<Scalar>
43 create_DfDp_mv(
44  const ModelEvaluator<Scalar>& model,
45  int l,
47  );
48 
49 
51 template<class Scalar>
52 ModelEvaluatorBase::DerivativeMultiVector<Scalar>
53 create_DgDx_dot_mv(
54  const ModelEvaluator<Scalar>& model,
55  int j,
57  );
58 
59 
61 template<class Scalar>
62 ModelEvaluatorBase::DerivativeMultiVector<Scalar>
63 create_DgDx_mv(
64  const ModelEvaluator<Scalar>& model,
65  int j,
67  );
68 
69 
71 template<class Scalar>
72 ModelEvaluatorBase::DerivativeMultiVector<Scalar>
73 create_DgDp_mv(
74  const ModelEvaluator<Scalar>& model,
75  int j,
76  int l,
78  );
79 
80 
82 template<class Scalar>
83 ModelEvaluatorBase::DerivativeMultiVector<Scalar>
84 get_dmv(
85  const ModelEvaluatorBase::Derivative<Scalar> &deriv
86  ,const std::string &derivName
87  );
88 
89 
91 template<class Scalar>
92 RCP<MultiVectorBase<Scalar> >
93 get_mv(
94  const ModelEvaluatorBase::Derivative<Scalar> &deriv
95  ,const std::string &derivName
97  );
98 
99 
105 template<class Scalar>
106 void assertDerivSpaces(
107  const std::string &modelEvalDescription,
108  const ModelEvaluatorBase::Derivative<Scalar> &deriv,
109  const std::string &deriv_name,
110  const VectorSpaceBase<Scalar> &fnc_space,
111  const std::string &fnc_space_name,
112  const VectorSpaceBase<Scalar> &var_space,
113  const std::string &var_space_name
114  );
115 
116 
121 template<class Scalar>
122 void assertInArgsOutArgsSetup(
123  const std::string &modelEvalDescription,
124  const ModelEvaluatorBase::InArgs<Scalar> &inArgs,
125  const ModelEvaluatorBase::OutArgs<Scalar> &outArgs
126  );
127 
128 
133 template<class Scalar>
134 void assertInArgsEvalObjects(
135  const ModelEvaluator<Scalar> &model,
136  const ModelEvaluatorBase::InArgs<Scalar> &inArgs
137  );
138 
139 
144 template<class Scalar>
145 void assertOutArgsEvalObjects(
146  const ModelEvaluator<Scalar> &model,
147  const ModelEvaluatorBase::OutArgs<Scalar> &outArgs,
148  const ModelEvaluatorBase::InArgs<Scalar> *inArgs = 0
149  );
150 
151 
153 template<class Scalar>
154 void eval_f(
155  const ModelEvaluator<Scalar> &model
156  ,const VectorBase<Scalar> &x
157  ,VectorBase<Scalar> *f
158  );
159 
160 
162 template<class Scalar>
163 void eval_f_W(
164  const ModelEvaluator<Scalar> &model
165  ,const VectorBase<Scalar> &x
166  ,VectorBase<Scalar> *f
167  ,LinearOpWithSolveBase<Scalar> *W
168  );
169 
170 
172 template<class Scalar>
173 void eval_f(
174  const ModelEvaluator<Scalar> &model
175  ,const VectorBase<Scalar> &x
176  ,const Scalar &t
177  ,VectorBase<Scalar> *f
178  );
179 
180 
182 template<class Scalar>
183 void eval_g(
184  const ModelEvaluator<Scalar> &model,
185  const int l,
186  const VectorBase<Scalar> &p_l,
187  const int j,
188  const Ptr<VectorBase<Scalar> > &g_j
189  );
190 
191 
193 template<class Scalar>
194 void eval_g(
195  const ModelEvaluator<Scalar> &model,
196  const int l,
197  const VectorBase<Scalar> &p_l,
198  const Scalar &t,
199  const int j,
200  VectorBase<Scalar> *g_j
201  );
202 
203 
205 template<class Scalar>
206 void eval_g_DgDp(
207  const ModelEvaluator<Scalar> &model,
208  const int l,
209  const VectorBase<Scalar> &p_l,
210  const int j,
211  const Ptr<VectorBase<Scalar> > &g_j,
212  const ModelEvaluatorBase::Derivative<Scalar> &DgDp_j_l
213  );
214 
215 
217 template<class Scalar>
218 void eval_f(
219  const ModelEvaluator<Scalar> &model
220  ,const VectorBase<Scalar> &x_dot
221  ,const VectorBase<Scalar> &x
222  ,const typename ModelEvaluatorBase::InArgs<Scalar>::ScalarMag &t
223  ,VectorBase<Scalar> *f
224  );
225 
226 
229 template<class Scalar>
230 void eval_f_W(
231  const ModelEvaluator<Scalar> &model
232  ,const VectorBase<Scalar> &x_dot
233  ,const VectorBase<Scalar> &x
234  ,const typename ModelEvaluatorBase::InArgs<Scalar>::ScalarMag &t
235  ,const Scalar &alpha
236  ,const Scalar &beta
237  ,VectorBase<Scalar> *f
238  ,LinearOpWithSolveBase<Scalar> *W
239  );
240 
241 
242 #ifdef HAVE_THYRA_ME_POLYNOMIAL
243 
244 
246 template<class Scalar>
247 void eval_f_poly(
248  const ModelEvaluator<Scalar> &model
249  ,const Teuchos::Polynomial< VectorBase<Scalar> > &x_poly
250  ,const typename ModelEvaluatorBase::InArgs<Scalar>::ScalarMag &t
251  ,Teuchos::Polynomial< VectorBase<Scalar> > *f_poly
252  );
253 
254 
256 template<class Scalar>
257 void eval_f_poly(
258  const ModelEvaluator<Scalar> &model
259  ,const Teuchos::Polynomial< VectorBase<Scalar> > &x_dot_poly
260  ,const VectorBase<Scalar> &x_poly
261  ,const typename ModelEvaluatorBase::InArgs<Scalar>::ScalarMag &t
262  ,Teuchos::Polynomial< VectorBase<Scalar> > *f_poly
263  );
264 
265 
266 #endif // HAVE_THYRA_ME_POLYNOMIAL
267 
268 
269 } // namespace Thyra
270 
271 
272 //
273 // Implementations
274 //
275 
276 
277 #include "Thyra_AssertOp.hpp"
278 #include "Teuchos_Utils.hpp"
279 
280 
281 template<class Scalar>
283 Thyra::clone( const ModelEvaluatorBase::InArgs<Scalar> &inArgs )
284 {
285  RCP<ModelEvaluatorBase::InArgs<Scalar> >
286  newInArgs = Teuchos::rcp(new ModelEvaluatorBase::InArgs<Scalar>);
287  *newInArgs = inArgs;
288  return newInArgs;
289 }
290 
291 
292 template<class Scalar>
294 Thyra::derivativeGradient(
295  const RCP<MultiVectorBase<Scalar> > &grad
296  )
297 {
298  return ModelEvaluatorBase::Derivative<Scalar>(
299  grad,
300  ModelEvaluatorBase::DERIV_MV_GRADIENT_FORM
301  );
302 }
303 
304 
305 template<class Scalar>
307 Thyra::create_DfDp_mv(
308  const ModelEvaluator<Scalar>& model,
309  int l,
310  ModelEvaluatorBase::EDerivativeMultiVectorOrientation orientation
311  )
312 {
313  TEUCHOS_TEST_FOR_EXCEPT(!(orientation==ModelEvaluatorBase::DERIV_MV_BY_COL));
314  return createMembers( model.get_f_space(), model.get_p_space(l)->dim() );
315 }
316 
317 
318 template<class Scalar>
320 Thyra::create_DgDx_dot_mv(
321  const ModelEvaluator<Scalar>& model,
322  int j,
323  ModelEvaluatorBase::EDerivativeMultiVectorOrientation orientation
324  )
325 {
326  typedef ModelEvaluatorBase MEB;
327  switch(orientation) {
328  case MEB::DERIV_MV_BY_COL:
329  return
330  MEB::DerivativeMultiVector<Scalar>(
331  createMembers( model.get_g_space(j), model.get_x_space()->dim() )
332  ,MEB::DERIV_MV_BY_COL
333  );
334  case MEB::DERIV_TRANS_MV_BY_ROW:
335  return
336  MEB::DerivativeMultiVector<Scalar>(
337  createMembers( model.get_x_space(), model.get_g_space(j)->dim() )
338  ,MEB::DERIV_TRANS_MV_BY_ROW
339  );
340  default:
342  }
343  TEUCHOS_UNREACHABLE_RETURN(MEB::DerivativeMultiVector<Scalar>());
344 }
345 
346 
347 template<class Scalar>
349 Thyra::create_DgDx_mv(
350  const ModelEvaluator<Scalar>& model,
351  int j,
352  ModelEvaluatorBase::EDerivativeMultiVectorOrientation orientation
353  )
354 {
355  return create_DgDx_dot_mv(model,j,orientation);
356 }
357 
358 
359 template<class Scalar>
361 Thyra::create_DgDp_mv(
362  const ModelEvaluator<Scalar>& model,
363  int j,
364  int l,
365  ModelEvaluatorBase::EDerivativeMultiVectorOrientation orientation
366  )
367 {
368  typedef ModelEvaluatorBase MEB;
369  switch(orientation) {
370  case MEB::DERIV_MV_BY_COL:
371  return
372  MEB::DerivativeMultiVector<Scalar>(
373  createMembers( model.get_g_space(j), model.get_p_space(l)->dim() )
374  ,MEB::DERIV_MV_BY_COL
375  );
376  case MEB::DERIV_TRANS_MV_BY_ROW:
377  return
378  MEB::DerivativeMultiVector<Scalar>(
379  createMembers( model.get_p_space(l), model.get_g_space(j)->dim() )
380  ,MEB::DERIV_TRANS_MV_BY_ROW
381  );
382  default:
384  }
385  TEUCHOS_UNREACHABLE_RETURN(MEB::DerivativeMultiVector<Scalar>());
386 }
387 
388 
389 template<class Scalar>
391 Thyra::get_dmv(
392  const ModelEvaluatorBase::Derivative<Scalar> &deriv
393  ,const std::string &derivName
394  )
395 {
397  deriv.getLinearOp().get()!=NULL, std::logic_error
398  ,"Error, LinearOpBase type not expected for " << derivName <<"!"
399  );
400  return deriv.getDerivativeMultiVector();
401 }
402 
403 
404 template<class Scalar>
406 Thyra::get_mv(
407  const ModelEvaluatorBase::Derivative<Scalar> &deriv
408  ,const std::string &derivName
409  ,ModelEvaluatorBase::EDerivativeMultiVectorOrientation orientation
410  )
411 {
412  typedef ModelEvaluatorBase MEB;
414  deriv.getLinearOp().get()!=NULL, std::logic_error
415  ,"Error, LinearOpBase type not expected for " << derivName <<"!"
416  );
417  MEB::DerivativeMultiVector<Scalar>
418  dmv = deriv.getDerivativeMultiVector();
419  RCP<MultiVectorBase<Scalar> >
420  mv = dmv.getMultiVector();
421  if( mv.get() ) {
423  dmv.getOrientation() != orientation, std::logic_error
424  ,"Error, the orientation " << toString(dmv.getOrientation()) << " is not the"
425  " expected orientation of " << toString(orientation)
426  << " for " << derivName << "!"
427  );
428  }
429  return mv;
430 }
431 
432 
433 template<class Scalar>
434 void Thyra::assertDerivSpaces(
435  const std::string &modelEvalDescription,
436  const ModelEvaluatorBase::Derivative<Scalar> &deriv,
437  const std::string &deriv_name,
438  const VectorSpaceBase<Scalar> &fnc_space,
439  const std::string &fnc_space_name,
440  const VectorSpaceBase<Scalar> &var_space,
441  const std::string &var_space_name
442  )
443 {
444  typedef ModelEvaluatorBase MEB;
445  if (!is_null(deriv.getLinearOp())) {
446  const RCP<const LinearOpBase<Scalar> > lo = deriv.getLinearOp();
447  if (!is_null(lo->range())) {
449  modelEvalDescription,
450  *lo->range(), deriv_name + ".range()",
451  fnc_space, fnc_space_name
452  );
454  modelEvalDescription,
455  *lo->domain(), deriv_name + ".domain()",
456  var_space, var_space_name
457  );
458  }
459  }
460  else if(!is_null(deriv.getMultiVector())) {
461  const RCP<const LinearOpBase<Scalar> > mv = deriv.getMultiVector();
462  switch(deriv.getMultiVectorOrientation()) {
463  case MEB::DERIV_MV_BY_COL: {
465  modelEvalDescription,
466  *mv->range(), deriv_name + ".range()",
467  fnc_space, fnc_space_name
468  );
470  modelEvalDescription,
471  *mv->domain(), deriv_name + ".domain()",
472  var_space, var_space_name
473  );
474  break;
475  }
476  case MEB::DERIV_TRANS_MV_BY_ROW: {
478  modelEvalDescription,
479  *mv->range(), deriv_name + "^T.range()",
480  var_space, var_space_name
481  );
483  modelEvalDescription,
484  *mv->domain(), deriv_name + "^T.domain()",
485  fnc_space, fnc_space_name
486  );
487  break;
488  }
489 #ifdef TEUCHOS_DEBUG
490  default:
492 #endif
493  }
494  }
495 }
496 
497 
498 template<class Scalar>
499 void Thyra::assertInArgsOutArgsSetup(
500  const std::string &modelEvalDescription,
501  const ModelEvaluatorBase::InArgs<Scalar> &inArgs,
502  const ModelEvaluatorBase::OutArgs<Scalar> &outArgs
503  )
504 {
505 
506  typedef ModelEvaluatorBase MEB;
507 
508  const int Ng = outArgs.Ng();
509  const int Np = outArgs.Np();
510 
511  // Description
512  TEUCHOS_ASSERT_EQUALITY(inArgs.modelEvalDescription(), modelEvalDescription);
513  TEUCHOS_ASSERT_EQUALITY(outArgs.modelEvalDescription(), modelEvalDescription);
514 
515  // Np
517  inArgs.Np() != outArgs.Np(), std::logic_error,
518  "Error: The underlying model " << modelEvalDescription << " incorrectly\n"
519  "set inArgs.Np() = "<<inArgs.Np()<<" != outArgs.Np() = "
520  <<outArgs.Np()<<"!"
521  );
522 
523  // x_dot
525  inArgs.supports(MEB::IN_ARG_x_dot) && !inArgs.supports(MEB::IN_ARG_x),
526  std::logic_error,
527  "Error: The underlying model " << modelEvalDescription << " supports\n"
528  "x_dot but does not support x!"
529  );
530 
531  // t
533  inArgs.supports(MEB::IN_ARG_x_dot) && !inArgs.supports(MEB::IN_ARG_t),
534  std::logic_error,
535  "Error: The underlying model " << modelEvalDescription << " supports\n"
536  "x_dot but does not support t!"
537  );
538 
539  // W and W_op
541  (
542  ( outArgs.supports(MEB::OUT_ARG_W) || outArgs.supports(MEB::OUT_ARG_W_op) )
543  &&
544  !inArgs.supports(MEB::IN_ARG_x)
545  ),
546  std::logic_error,
547  "Error: The underlying model " << modelEvalDescription << " says that\n"
548  "it supports W and/or W_op but it does not support x!"
549  );
551  (
552  ( outArgs.supports(MEB::OUT_ARG_W) || outArgs.supports(MEB::OUT_ARG_W_op) )
553  &&
554  inArgs.supports(MEB::IN_ARG_x_dot)
555  &&
556  !( inArgs.supports(MEB::IN_ARG_alpha) && inArgs.supports(MEB::IN_ARG_beta) )
557  ),
558  std::logic_error,
559  "Error: The underlying model " << modelEvalDescription << " supports W and/or W_op\n"
560  "and x_dot but it does not support alpha and beta as InArgs! \n"
561  "If the model supports W and x_dot, then it can be interpreted as an implicit \n"
562  "ODE/DAE and therefore D(f)/D(x_dot) must be non-zero and therefore alpha must \n"
563  "be supported. If, however, the model can be interpreted as an explicit ODE \n"
564  "then x_dot should not be supported at all."
565  );
566 
567  for ( int l = 0; l < Np; ++l ) {
568 
569  // DfDp(l): OutArgs checks this automatically!
570 
571  for ( int j = 0; j < Ng; ++j ) {
572 
573  // DgDx_dot(j)
575  ( !outArgs.supports(MEB::OUT_ARG_DgDx_dot,j).none()
576  && !inArgs.supports(MEB::IN_ARG_x_dot) ),
577  std::logic_error,
578  "Error: The underlying model " << modelEvalDescription << " says that\n"
579  "it supports DgDx_dot("<<j<<") but it does not support x_dot!"
580  );
581 
582  // DgDx(j)
584  ( !outArgs.supports(MEB::OUT_ARG_DgDx,j).none()
585  && !inArgs.supports(MEB::IN_ARG_x) ),
586  std::logic_error,
587  "Error: The underlying model " << modelEvalDescription << " says that\n"
588  "it supports DgDx("<<j<<") but it does not support x!"
589  );
590 
591  // DgDp(j,l): OutArgs checks this automatically!
592 
593  }
594 
595  }
596 
597 }
598 
599 
600 template<class Scalar>
601 void Thyra::assertInArgsEvalObjects(
602  const ModelEvaluator<Scalar> &model,
603  const ModelEvaluatorBase::InArgs<Scalar> &inArgs
604  )
605 {
606 
607  typedef ModelEvaluatorBase MEB;
608 
609  const std::string description = model.description();
610  const int Np = inArgs.Np();
611 
612  model.createInArgs().assertSameSupport(inArgs);
613 
614  // x_dot
615  if ( inArgs.supports(MEB::IN_ARG_x_dot) && !is_null(inArgs.get_x_dot()) ) {
617  description, *inArgs.get_x_dot()->space(), *model.get_x_space() );
618  }
619 
620  // x
621  if ( inArgs.supports(MEB::IN_ARG_x) && !is_null(inArgs.get_x()) ) {
623  description, *inArgs.get_x()->space(), *model.get_x_space() );
624  }
625 
626  // p(l)
627  for ( int l = 0; l < Np; ++l ) {
628  if (!is_null(inArgs.get_p(l))) {
630  description, *inArgs.get_p(l)->space(), *model.get_p_space(l) );
631  }
632  }
633 
634 }
635 
636 
637 template<class Scalar>
638 void Thyra::assertOutArgsEvalObjects(
639  const ModelEvaluator<Scalar> &model,
640  const ModelEvaluatorBase::OutArgs<Scalar> &outArgs,
641  const ModelEvaluatorBase::InArgs<Scalar> *inArgs
642  )
643 {
644 
645  typedef ScalarTraits<Scalar> ST;
646  typedef Teuchos::Utils TU;
647  typedef ModelEvaluatorBase MEB;
648 
649  const std::string description = model.description();
650  const int Ng = outArgs.Ng();
651  const int Np = outArgs.Np();
652 
653  if (inArgs) {
654  TEUCHOS_ASSERT_EQUALITY(outArgs.Np(), inArgs->Np());
655  }
656 
657  model.createOutArgs().assertSameSupport(outArgs);
658 
659  // f
660  if ( outArgs.supports(MEB::OUT_ARG_f) && !is_null(outArgs.get_f()) ) {
662  description, *outArgs.get_f()->space(), *model.get_f_space() );
663  }
664 
665  // W
666  if ( outArgs.supports(MEB::OUT_ARG_W) && !is_null(outArgs.get_W()) ) {
667  if (!is_null(outArgs.get_W()->range())) {
669  description, *outArgs.get_W()->range(), *model.get_f_space() );
671  description, *outArgs.get_W()->domain(), *model.get_x_space() );
672  }
673  }
674 
675  // W_op
676  if ( outArgs.supports(MEB::OUT_ARG_W_op) && !is_null(outArgs.get_W_op()) ) {
677  if (!is_null(outArgs.get_W_op()->range())) {
679  description, *outArgs.get_W_op()->range(), *model.get_f_space() );
681  description, *outArgs.get_W_op()->domain(), *model.get_x_space() );
682  }
683  }
684 
685  // alpha and beta (not really in outArgs but can only be validated if W or
686  // W_op is set)
687  if (
688  inArgs
689  &&
690  (
691  ( outArgs.supports(MEB::OUT_ARG_W) && !is_null(outArgs.get_W()) )
692  ||
693  ( outArgs.supports(MEB::OUT_ARG_W_op) && !is_null(outArgs.get_W_op()) )
694  )
695  )
696  {
697  if ( inArgs->supports(MEB::IN_ARG_alpha) && inArgs->supports(MEB::IN_ARG_beta) ) {
698  // 08/25/08 tscoffe: In the block-composed linear operator case for
699  // Rythmos::ImplicitRKModelEvaluator, I need to specify that a given
700  // block is all zeros and I'm depending on the underlying model to
701  // intelligently fill the block with zeros if both alpha and beta are
702  // zero.
703  //TEUCHOS_TEST_FOR_EXCEPT( inArgs->get_alpha() == ST::zero() && inArgs->get_beta() == ST::zero() );
704  }
705  else if ( inArgs->supports(MEB::IN_ARG_beta) ) {
706  TEUCHOS_TEST_FOR_EXCEPT( inArgs->get_beta() == ST::zero() );
707  }
708  }
709 
710  // DfDp(l)
711  if (outArgs.supports(MEB::OUT_ARG_f)) {
712  for ( int l = 0; l < Np; ++l ) {
713  if (!outArgs.supports(MEB::OUT_ARG_DfDp,l).none()) {
715  description,
716  outArgs.get_DfDp(l), "DfDp("+TU::toString(l)+")",
717  *model.get_f_space(), "f_space",
718  *model.get_p_space(l), "p_space("+TU::toString(l)+")"
719  );
720  }
721  }
722  }
723 
724  // g(l)
725  for ( int j = 0; j < Ng; ++j ) {
726  if (!is_null(outArgs.get_g(j))) {
728  description, *outArgs.get_g(j)->space(), *model.get_g_space(j) );
729  }
730  }
731 
732  // DgDx_dot(j)
733  for ( int j = 0; j < Ng; ++j ) {
734  if (!outArgs.supports(MEB::OUT_ARG_DgDx_dot,j).none()) {
736  description,
737  outArgs.get_DgDx_dot(j), "DgDx_dot("+TU::toString(j)+")",
738  *model.get_g_space(j), "g_space("+TU::toString(j)+")",
739  *model.get_x_space(), "x_space"
740  );
741  }
742  }
743 
744  // DgDx(j)
745  for ( int j = 0; j < Ng; ++j ) {
746  if (!outArgs.supports(MEB::OUT_ARG_DgDx,j).none()) {
748  description,
749  outArgs.get_DgDx(j), "DgDx("+TU::toString(j)+")",
750  *model.get_g_space(j), "g_space("+TU::toString(j)+")",
751  *model.get_x_space(), "x_space"
752  );
753  }
754  }
755 
756  // Assert DgDp(j,l)
757  for ( int j = 0; j < Ng; ++j ) {
758  for ( int l = 0; l < Np; ++l ) {
759  if (!outArgs.supports(MEB::OUT_ARG_DgDp,j,l).none()) {
760  const std::string j_str = TU::toString(j);
761  const std::string l_str = TU::toString(l);
763  description,
764  outArgs.get_DgDp(j,l), "DgDp("+j_str+","+l_str+")",
765  *model.get_g_space(j), "g_space("+j_str+")",
766  *model.get_p_space(l), "p_space("+l_str+")"
767  );
768  }
769  }
770  }
771 
772 }
773 
774 
775 template<class Scalar>
776 void Thyra::eval_f(
777  const ModelEvaluator<Scalar> &model
778  ,const VectorBase<Scalar> &x
779  ,VectorBase<Scalar> *f
780  )
781 {
782  typedef ModelEvaluatorBase MEB;
783  MEB::InArgs<Scalar> inArgs = model.createInArgs();
784  MEB::OutArgs<Scalar> outArgs = model.createOutArgs();
785  inArgs.set_x(Teuchos::rcp(&x,false));
786  outArgs.set_f(Teuchos::rcp(f,false));
787  model.evalModel(inArgs,outArgs);
788 }
789 
790 
791 template<class Scalar>
792 void Thyra::eval_f_W(
793  const ModelEvaluator<Scalar> &model
794  ,const VectorBase<Scalar> &x
795  ,VectorBase<Scalar> *f
796  ,LinearOpWithSolveBase<Scalar> *W
797  )
798 {
799 
800  typedef ModelEvaluatorBase MEB;
801 
802  MEB::InArgs<Scalar> inArgs = model.createInArgs();
803  MEB::OutArgs<Scalar> outArgs = model.createOutArgs();
804 
805  inArgs.set_x(Teuchos::rcp(&x,false));
806 
807  if (f) outArgs.set_f(Teuchos::rcp(f,false));
808  if (W) outArgs.set_W(Teuchos::rcp(W,false));
809 
810  model.evalModel(inArgs,outArgs);
811 
812 }
813 
814 
815 template<class Scalar>
816 void Thyra::eval_f(
817  const ModelEvaluator<Scalar> &model
818  ,const VectorBase<Scalar> &x
819  ,const Scalar &t
820  ,VectorBase<Scalar> *f
821  )
822 {
823  typedef ModelEvaluatorBase MEB;
824  MEB::InArgs<Scalar> inArgs = model.createInArgs();
825  MEB::OutArgs<Scalar> outArgs = model.createOutArgs();
826  inArgs.set_x(Teuchos::rcp(&x,false));
827  if(inArgs.supports(MEB::IN_ARG_t)) inArgs.set_t(t);
828  outArgs.set_f(Teuchos::rcp(f,false));
829  model.evalModel(inArgs,outArgs);
830 }
831 
832 
833 template<class Scalar>
834 void Thyra::eval_g(
835  const ModelEvaluator<Scalar> &model,
836  const int l,
837  const VectorBase<Scalar> &p_l,
838  const int j,
839  const Ptr<VectorBase<Scalar> > &g_j
840  )
841 {
842  typedef ModelEvaluatorBase MEB;
843  MEB::InArgs<Scalar> inArgs = model.createInArgs();
844  MEB::OutArgs<Scalar> outArgs= model.createOutArgs();
845  inArgs.set_p(l, Teuchos::rcpFromRef(p_l));
846  outArgs.set_g(j, Teuchos::rcpFromRef(*g_j));
847  model.evalModel(inArgs,outArgs);
848 }
849 
850 
851 template<class Scalar>
852 void Thyra::eval_g(
853  const ModelEvaluator<Scalar> &model,
854  const int l,
855  const VectorBase<Scalar> &p_l,
856  const Scalar &t,
857  const int j,
858  VectorBase<Scalar> *g_j
859  )
860 {
861  typedef ModelEvaluatorBase MEB;
862  MEB::InArgs<Scalar> inArgs = model.createInArgs();
863  MEB::OutArgs<Scalar> outArgs= model.createOutArgs();
864  inArgs.set_p(l,Teuchos::rcp(&p_l,false));
865  inArgs.set_t(t);
866  outArgs.set_g(j,Teuchos::rcp(g_j,false));
867  model.evalModel(inArgs,outArgs);
868 }
869 
870 
871 template<class Scalar>
872 void Thyra::eval_g_DgDp(
873  const ModelEvaluator<Scalar> &model,
874  const int l,
875  const VectorBase<Scalar> &p_l,
876  const int j,
877  const Ptr<VectorBase<Scalar> > &g_j,
878  const ModelEvaluatorBase::Derivative<Scalar> &DgDp_j_l
879  )
880 {
881  typedef ModelEvaluatorBase MEB;
882  MEB::InArgs<Scalar> inArgs = model.createInArgs();
883  MEB::OutArgs<Scalar> outArgs= model.createOutArgs();
884  inArgs.set_p(l, Teuchos::rcpFromRef(p_l));
885  if (!is_null(g_j)) {
886  outArgs.set_g(j, Teuchos::rcpFromPtr(g_j));
887  }
888  if (!DgDp_j_l.isEmpty()) {
889  outArgs.set_DgDp(j, l, DgDp_j_l);
890  }
891  model.evalModel(inArgs,outArgs);
892 }
893 
894 
895 template<class Scalar>
896 void Thyra::eval_f(
897  const ModelEvaluator<Scalar> &model
898  ,const VectorBase<Scalar> &x_dot
899  ,const VectorBase<Scalar> &x
900  ,const typename ModelEvaluatorBase::InArgs<Scalar>::ScalarMag &t
901  ,VectorBase<Scalar> *f
902  )
903 {
904 
905  typedef ModelEvaluatorBase MEB;
906 
907  MEB::InArgs<Scalar> inArgs = model.createInArgs();
908  MEB::OutArgs<Scalar> outArgs = model.createOutArgs();
909 
910  inArgs.set_x_dot(Teuchos::rcp(&x_dot,false));
911  inArgs.set_x(Teuchos::rcp(&x,false));
912  if(inArgs.supports(MEB::IN_ARG_t))
913  inArgs.set_t(t);
914 
915  outArgs.set_f(Teuchos::rcp(f,false));
916 
917  model.evalModel(inArgs,outArgs);
918 
919 }
920 
921 
922 template<class Scalar>
923 void Thyra::eval_f_W(
924  const ModelEvaluator<Scalar> &model
925  ,const VectorBase<Scalar> &x_dot
926  ,const VectorBase<Scalar> &x
927  ,const typename ModelEvaluatorBase::InArgs<Scalar>::ScalarMag &t
928  ,const Scalar &alpha
929  ,const Scalar &beta
930  ,VectorBase<Scalar> *f
931  ,LinearOpWithSolveBase<Scalar> *W
932  )
933 {
934 
935  typedef ModelEvaluatorBase MEB;
936 
937  MEB::InArgs<Scalar> inArgs = model.createInArgs();
938  MEB::OutArgs<Scalar> outArgs = model.createOutArgs();
939 
940  inArgs.set_x_dot(Teuchos::rcp(&x_dot,false));
941  inArgs.set_x(Teuchos::rcp(&x,false));
942  if(inArgs.supports(MEB::IN_ARG_t))
943  inArgs.set_t(t);
944  inArgs.set_alpha(alpha);
945  inArgs.set_beta(beta);
946 
947  if(f) outArgs.set_f(Teuchos::rcp(f,false));
948  if(W) outArgs.set_W(Teuchos::rcp(W,false));
949 
950  model.evalModel(inArgs,outArgs);
951 
952 }
953 
954 
955 #ifdef HAVE_THYRA_ME_POLYNOMIAL
956 
957 
958 template<class Scalar>
959 void Thyra::eval_f_poly(
960  const ModelEvaluator<Scalar> &model
961  ,const Teuchos::Polynomial< VectorBase<Scalar> > &x_poly
962  ,const typename ModelEvaluatorBase::InArgs<Scalar>::ScalarMag &t
963  ,Teuchos::Polynomial< VectorBase<Scalar> > *f_poly
964  )
965 {
966 
967  typedef ModelEvaluatorBase MEB;
968 
969  MEB::InArgs<Scalar> inArgs = model.createInArgs();
970  MEB::OutArgs<Scalar> outArgs = model.createOutArgs();
971 
972  inArgs.set_x_poly(Teuchos::rcp(&x_poly,false));
973  if(inArgs.supports(MEB::IN_ARG_t))
974  inArgs.set_t(t);
975 
976  outArgs.set_f_poly(Teuchos::rcp(f_poly,false));
977 
978  model.evalModel(inArgs,outArgs);
979 
980 }
981 
982 
983 template<class Scalar>
984 void Thyra::eval_f_poly(
985  const ModelEvaluator<Scalar> &model
986  ,const Teuchos::Polynomial< VectorBase<Scalar> > &x_dot_poly
987  ,const VectorBase<Scalar> &x_poly
988  ,const typename ModelEvaluatorBase::InArgs<Scalar>::ScalarMag &t
989  ,Teuchos::Polynomial< VectorBase<Scalar> > *f_poly
990  )
991 {
992 
993  typedef ModelEvaluatorBase MEB;
994 
995  MEB::InArgs<Scalar> inArgs = model.createInArgs();
996  MEB::OutArgs<Scalar> outArgs = model.createOutArgs();
997 
998  inArgs.set_x_dot_poly(Teuchos::rcp(&x_dot_poly,false));
999  inArgs.set_x_poly(Teuchos::rcp(&x_poly,false));
1000  if(inArgs.supports(MEB::IN_ARG_t))
1001  inArgs.set_t(t);
1002 
1003  outArgs.set_f_poly(Teuchos::rcp(f_poly,false));
1004 
1005  model.evalModel(inArgs,outArgs);
1006 
1007 }
1008 
1009 
1010 #endif // HAVE_THYRA_ME_POLYNOMIAL
1011 
1012 
1013 #endif // THYRA_MODEL_EVALUATOR_HELPERS_HPP
#define THYRA_ASSERT_VEC_SPACES(FUNC_NAME, VS1, VS2)
This is a very useful macro that should be used to validate that two vector spaces are compatible...
bool is_null(const boost::shared_ptr< T > &p)
std::string toString(ModelEvaluatorBase::EInArgsMembers)
#define TEUCHOS_TEST_FOR_EXCEPTION(throw_exception_test, Exception, msg)
Simple aggregate class that stores a derivative object as a general linear operator or as a multi-vec...
Simple aggregate class for a derivative object represented as a column-wise multi-vector or its trans...
TEUCHOS_DEPRECATED RCP< T > rcp(T *p, Dealloc_T dealloc, bool owns_mem)
TEUCHOSCORE_LIB_DLL_EXPORT std::string toString(const EVerbosityLevel verbLevel)
ModelEvaluatorBase::DerivativeMultiVector< Scalar > create_DgDx_dot_mv(const ModelEvaluator< Scalar > &model, int j, ModelEvaluatorBase::EDerivativeMultiVectorOrientation orientation)
#define TEUCHOS_UNREACHABLE_RETURN(dummyReturnVal)
void assertDerivSpaces(const std::string &modelEvalDescription, const ModelEvaluatorBase::Derivative< Scalar > &deriv, const std::string &deriv_name, const VectorSpaceBase< Scalar > &fnc_space, const std::string &fnc_space_name, const VectorSpaceBase< Scalar > &var_space, const std::string &var_space_name)
Assert that that Thyra objects imbedded in a Derivative object matches its function and variable spac...
#define THYRA_ASSERT_VEC_SPACES_NAMES(FUNC_NAME, VS1, VS1_NAME, VS2, VS2_NAME)
Helper assertion macro.
#define TEUCHOS_ASSERT_EQUALITY(val1, val2)
RCP< MultiVectorBase< Scalar > > createMembers(const RCP< const VectorSpaceBase< Scalar > > &vs, int numMembers, const std::string &label="")
Create a set of vector members (a MultiVectorBase) from the vector space.
#define TEUCHOS_TEST_FOR_EXCEPT(throw_exception_test)