Thyra  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Thyra_DefaultFiniteDifferenceModelEvaluator_def.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_DEFAULT_FINITE_DIFFERENCE_MODEL_EVALUATOR_DEF_HPP
11 #define THYRA_DEFAULT_FINITE_DIFFERENCE_MODEL_EVALUATOR_DEF_HPP
12 
13 #include "Thyra_DefaultFiniteDifferenceModelEvaluator_decl.hpp"
14 
15 
16 namespace Thyra {
17 
18 
19 // Constructors/initializers/accessors/utilities
20 
21 
22 template<class Scalar>
24 {}
25 
26 
27 template<class Scalar>
29  const RCP<ModelEvaluator<Scalar> > &thyraModel,
30  const RCP<DirectionalFiniteDiffCalculator<Scalar> > &direcFiniteDiffCalculator_in
31  )
32 {
34  direcFiniteDiffCalculator_ = direcFiniteDiffCalculator_in;
35 }
36 
37 
38 // Public functions overridden from Teuchos::Describable
39 
40 
41 template<class Scalar>
43 {
45  thyraModel = this->getUnderlyingModel();
46  std::ostringstream oss;
47  oss << "Thyra::DefaultFiniteDifferenceModelEvaluator{";
48  oss << "thyraModel=";
49  if(thyraModel.get())
50  oss << "\'"<<thyraModel->description()<<"\'";
51  else
52  oss << "NULL";
53  oss << "}";
54  return oss.str();
55 }
56 
57 
58 // Private functions overridden from ModelEvaulatorDefaultBase
59 
60 
61 template<class Scalar>
64 {
65  typedef ModelEvaluatorBase MEB;
67  thyraModel = this->getUnderlyingModel();
68  const MEB::OutArgs<Scalar> wrappedOutArgs = thyraModel->createOutArgs();
69  const int l_Np = wrappedOutArgs.Np(), l_Ng = wrappedOutArgs.Ng();
70  MEB::OutArgsSetup<Scalar> outArgs;
71  outArgs.setModelEvalDescription(this->description());
72  outArgs.set_Np_Ng(l_Np,l_Ng);
73  outArgs.setSupports(wrappedOutArgs);
74  if (wrappedOutArgs.supports(MEB::OUT_ARG_f)) {
75  for( int l = 0; l < l_Np; ++l ) {
76  outArgs.setSupports(MEB::OUT_ARG_DfDp,l,MEB::DERIV_MV_BY_COL);
77  }
78  }
79  for( int j = 0; j < l_Ng; ++j ) {
80  for( int l = 0; l < l_Np; ++l ) {
81  outArgs.setSupports( MEB::OUT_ARG_DgDp , j, l, MEB::DERIV_MV_BY_COL);
82  }
83  }
84  // ToDo: Add support for more derivatives as needed!
85  return outArgs;
86 }
87 
88 
89 template<class Scalar>
90 void DefaultFiniteDifferenceModelEvaluator<Scalar>::evalModelImpl(
91  const ModelEvaluatorBase::InArgs<Scalar> &inArgs,
92  const ModelEvaluatorBase::OutArgs<Scalar> &outArgs
93  ) const
94 {
95  using Teuchos::rcp;
96  using Teuchos::rcp_const_cast;
97  using Teuchos::rcp_dynamic_cast;
98  using Teuchos::OSTab;
99  typedef ModelEvaluatorBase MEB;
100  namespace DFDCT = DirectionalFiniteDiffCalculatorTypes;
101 
102  typedef RCP<VectorBase<Scalar> > V_ptr;
103 
104  THYRA_MODEL_EVALUATOR_DECORATOR_EVAL_MODEL_BEGIN(
105  "Thyra::DefaultFiniteDifferenceModelEvaluator",inArgs,outArgs
106  );
107 
108  //
109  // Note: Just do derivatives DfDp(l) and DgDp(j,l) for now!
110  //
111 
112  const RCP<const VectorSpaceBase<Scalar> >
113  p_space = thyraModel->get_p_space(0),
114  g_space = thyraModel->get_g_space(0);
115 
116  //
117  // A) Compute the base point
118  //
119 
120  if(out.get() && includesVerbLevel(verbLevel,Teuchos::VERB_LOW))
121  *out << "\nComputing the base point ...\n";
122 
123  const int l_Np = outArgs.Np();
124  const int l_Ng = outArgs.Ng();
125  MEB::InArgs<Scalar> wrappedInArgs = inArgs;
126  MEB::OutArgs<Scalar> baseFunc = thyraModel->createOutArgs();
127  if( outArgs.supports(MEB::OUT_ARG_f) && outArgs.get_f().get() )
128  baseFunc.set_f(outArgs.get_f());
129  for( int j = 0; j < l_Ng; ++j ) {
130  V_ptr g_j;
131  if( (g_j=outArgs.get_g(j)).get() )
132  baseFunc.set_g(j,g_j);
133  }
134  // 2007/08/27: We really should really try to allow some derivatives to pass
135  // through and some derivatives to be computed by finite differences. Right
136  // now, if you use this class, all derivatives w.r.t. parameters are finite
137  // differenced and that is not given the user enough control!
138 
139  thyraModel->evalModel(wrappedInArgs,baseFunc);
140 
141  bool failed = baseFunc.isFailed();
142 
143  //
144  // B) Compute the derivatives
145  //
146 
147  if(!failed) {
148 
149  // a) Determine what derivatives you need to support first
150 
151  Array<int> compute_DfDp;
152  Array<Array<int> > compute_DgDp(l_Ng);
153  DFDCT::SelectedDerivatives selectedDerivs;
154 
155  for ( int l = 0; l < l_Np; ++l ) {
156 
157  // DfDp(l)
158  if(
159  outArgs.supports(MEB::OUT_ARG_DfDp,l).none()==false
160  &&
161  outArgs.get_DfDp(l).isEmpty()==false
162  )
163  {
164  selectedDerivs.supports(MEB::OUT_ARG_DfDp,l);
165  compute_DfDp.push_back(true);
166  }
167  else
168  {
169  compute_DfDp.push_back(false);
170  }
171 
172  // DgDp(j=0...,l)
173  for ( int j = 0; j < l_Ng; ++j ) {
174  if(
175  outArgs.supports(MEB::OUT_ARG_DgDp,j,l).none()==false
176  &&
177  outArgs.get_DgDp(j,l).isEmpty()==false
178  )
179  {
180  selectedDerivs.supports(MEB::OUT_ARG_DgDp,j,l);
181  compute_DgDp[j].push_back(true);
182  }
183  else
184  {
185  compute_DgDp[j].push_back(false);
186  }
187  }
188  }
189 
190  // b) Create the deriv OutArgs and set the output objects that need to be
191  // computed with finite differences
192 
193  MEB::OutArgs<Scalar>
194  deriv = direcFiniteDiffCalculator_->createOutArgs(
195  *thyraModel, selectedDerivs );
196 
197  for ( int l = 0; l < l_Np; ++l ) {
198  if ( compute_DfDp[l] )
199  deriv.set_DfDp(l,outArgs.get_DfDp(l));
200  for ( int j = 0; j < l_Ng; ++j ) {
201  if ( compute_DgDp[j][l] )
202  deriv.set_DgDp(j,l,outArgs.get_DgDp(j,l));
203  }
204  }
205 
206  // c) Compute the missing functions with finite differences!
207 
208  direcFiniteDiffCalculator_->calcDerivatives(
209  *thyraModel,inArgs,baseFunc,deriv
210  );
211 
212  }
213 
214  if(failed) {
215  if(out.get() && includesVerbLevel(verbLevel,Teuchos::VERB_LOW))
216  *out << "\nEvaluation failed, returning NaNs ...\n";
217  outArgs.setFailed();
218  }
219 
220  THYRA_MODEL_EVALUATOR_DECORATOR_EVAL_MODEL_END();
221 
222 }
223 
224 
225 } // namespace Thyra
226 
227 
228 #endif // THYRA_DEFAULT_FINITE_DIFFERENCE_MODEL_EVALUATOR_DEF_HPP
Pure abstract base interface for evaluating a stateless &quot;model&quot; that can be mapped into a number of d...
basic_OSTab< char > OSTab
Concrete aggregate class for all output arguments computable by a ModelEvaluator subclass object...
void initialize(const RCP< ModelEvaluator< Scalar > > &model)
Initialize given a non-const model evaluator.
T * get() const
This class wraps any ModelEvaluator object and computes certain derivatives using finite differences...
TEUCHOS_DEPRECATED RCP< T > rcp(T *p, Dealloc_T dealloc, bool owns_mem)
TEUCHOSCORE_LIB_DLL_EXPORT bool includesVerbLevel(const EVerbosityLevel verbLevel, const EVerbosityLevel requestedVerbLevel, const bool isDefaultLevel=false)
Base subclass for ModelEvaluator that defines some basic types.
void initialize(const RCP< ModelEvaluator< Scalar > > &thyraModel, const RCP< DirectionalFiniteDiffCalculator< Scalar > > &direcFiniteDiffCalculator)
Utility class for computing directional finite differences of a model.