Thyra Version of the Day
Loading...
Searching...
No Matches
Thyra_EpetraModelEvaluator.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_EPETRA_MODEL_EVALUATOR_HPP
11#define THYRA_EPETRA_MODEL_EVALUATOR_HPP
12
13#include "Thyra_ModelEvaluatorDefaultBase.hpp"
15#include "Thyra_LinearOpWithSolveFactoryBase.hpp"
16#include "EpetraExt_ModelEvaluator.h"
17#include "Epetra_Map.h"
18#include "Teuchos_Array.hpp"
19
20
21namespace Thyra {
22
23
144 : public ModelEvaluatorDefaultBase<double>,
145 virtual public Teuchos::ParameterListAcceptor
146{
147public:
148
151
154
157 const RCP<const EpetraExt::ModelEvaluator> &epetraModel,
159 );
160
162 void initialize(
163 const RCP<const EpetraExt::ModelEvaluator> &epetraModel,
165 );
166
169
175 void setNominalValues( const ModelEvaluatorBase::InArgs<double>& nominalValues );
176
185 const RCP<const Epetra_Vector> &stateVariableScalingVec
186 );
187
192
197
201 const RCP<const Epetra_Vector> &stateFunctionScalingVec
202 );
203
208
210 void uninitialize(
211 RCP<const EpetraExt::ModelEvaluator> *epetraModel = NULL,
212 RCP<LinearOpWithSolveFactoryBase<double> > *W_factory = NULL
213 );
214
217
219 bool finalPointWasSolved() const;
220
222
225
227 std::string description() const;
228
230
233
235 void setParameterList(RCP<Teuchos::ParameterList> const& paramList);
244
246
249
251 int Np() const;
253 int Ng() const;
281 void reportFinalPoint(
282 const ModelEvaluatorBase::InArgs<double> &finalPoint
283 ,const bool wasSolved
284 );
285
287
288 // Made public to simplify implementation but this is harmless to be public.
289 // Clients should not deal with this type.
290 enum EStateFunctionScaling { STATE_FUNC_SCALING_NONE, STATE_FUNC_SCALING_ROW_SUM };
291
292private:
293
296
298 RCP<LinearOpBase<double> > create_DfDp_op_impl(int l) const;
300 RCP<LinearOpBase<double> > create_DgDx_dot_op_impl(int j) const;
302 RCP<LinearOpBase<double> > create_DgDx_op_impl(int j) const;
304 RCP<LinearOpBase<double> > create_DgDp_op_impl(int j, int l) const;
306 ModelEvaluatorBase::OutArgs<double> createOutArgsImpl() const;
308 void evalModelImpl(
311 ) const;
312
314
315private:
316
317 // ////////////////////
318 // Private types
319
322 typedef std::vector<bool> p_map_is_local_t;
323 typedef std::vector<bool> g_map_is_local_t;
324
326 p_space_t;
328 g_space_t;
329
330 // /////////////////////
331 // Private data members
332
334
336
338
340 p_map_t p_map_;
341 g_map_t g_map_;
342 p_map_is_local_t p_map_is_local_;
343 p_map_is_local_t g_map_is_local_;
345
347 p_space_t p_space_;
349 g_space_t g_space_;
350
351 mutable ModelEvaluatorBase::InArgs<double> nominalValues_;
352 mutable ModelEvaluatorBase::InArgs<double> lowerBounds_;
353 mutable ModelEvaluatorBase::InArgs<double> upperBounds_;
354 mutable bool nominalValuesAndBoundsAreUpdated_;
355
357
358 EStateFunctionScaling stateFunctionScaling_;
359 mutable RCP<const Epetra_Vector> stateFunctionScalingVec_;
360
361 RCP<const Epetra_Vector> stateVariableScalingVec_; // S_x
362 mutable RCP<const Epetra_Vector> invStateVariableScalingVec_; // inv(S_x)
363 mutable EpetraExt::ModelEvaluator::InArgs epetraInArgsScaling_;
364 mutable EpetraExt::ModelEvaluator::OutArgs epetraOutArgsScaling_;
365
366 mutable RCP<Epetra_Vector> x_unscaled_;
367 mutable RCP<Epetra_Vector> x_dot_unscaled_;
368
369 mutable ModelEvaluatorBase::InArgs<double> prototypeInArgs_;
370 mutable ModelEvaluatorBase::OutArgs<double> prototypeOutArgs_;
371 mutable bool currentInArgsOutArgs_;
372
373 bool finalPointWasSolved_;
374
375 // //////////////////////////
376 // Private member functions
377
379 void convertInArgsFromEpetraToThyra(
380 const EpetraExt::ModelEvaluator::InArgs &epetraInArgs,
382 ) const;
383
385 void convertInArgsFromThyraToEpetra(
387 EpetraExt::ModelEvaluator::InArgs *epetraInArgs
388 ) const;
389
391 void convertOutArgsFromThyraToEpetra(
392 // Thyra form of the outArgs
394 // Epetra form of the unscaled output arguments
395 EpetraExt::ModelEvaluator::OutArgs *epetraUnscaledOutArgs,
396 // The passed-in form of W
398 RCP<EpetraLinearOp> *efwdW,
399 // The actual Epetra object passed to the underylying EpetraExt::ModelEvaluator
401 ) const;
402
404 void preEvalScalingSetup(
405 EpetraExt::ModelEvaluator::InArgs *epetraInArgs,
406 EpetraExt::ModelEvaluator::OutArgs *epetraUnscaledOutArgs,
408 const Teuchos::EVerbosityLevel verbLevel
409 ) const;
410
412 void postEvalScalingSetup(
413 const EpetraExt::ModelEvaluator::OutArgs &epetraUnscaledOutArgs,
415 const Teuchos::EVerbosityLevel verbLevel
416 ) const;
417
419 void finishConvertingOutArgsFromEpetraToThyra(
420 const EpetraExt::ModelEvaluator::OutArgs &epetraOutArgs,
422 RCP<EpetraLinearOp> &efwdW,
424 const ModelEvaluatorBase::OutArgs<double> &outArgs // Output!
425 ) const;
426 // 2007/08/03: rabartl: Above, I pass many of the RCP objects by non-const
427 // reference since I don't want the compiler to perform any implicit
428 // conversions on this RCP objects.
429
431 void updateNominalValuesAndBounds() const;
432
434 void updateInArgsOutArgs() const;
435
437 RCP<EpetraLinearOp> create_epetra_W_op() const;
438
439};
440
441
442//
443// Utility functions
444//
445
446
452 const RCP<const EpetraExt::ModelEvaluator> &epetraModel,
454 );
455
456
461convert( const EpetraExt::ModelEvaluator::EDerivativeMultiVectorOrientation &mvOrientation );
462
463
467EpetraExt::ModelEvaluator::EDerivativeMultiVectorOrientation
469
470
475convert( const EpetraExt::ModelEvaluator::DerivativeProperties &derivativeProperties );
476
477
482convert( const EpetraExt::ModelEvaluator::DerivativeSupport &derivativeSupport );
483
484
488EpetraExt::ModelEvaluator::Derivative
491 const RCP<const Epetra_Map> &fnc_map,
492 const RCP<const Epetra_Map> &var_map
493 );
494
495EpetraExt::ModelEvaluator::MPDerivative
496convert(
497 const ModelEvaluatorBase::MPDerivative &derivative,
498 const RCP<const Epetra_Map> &fnc_map,
499 const RCP<const Epetra_Map> &var_map
500 );
501
502} // namespace Thyra
503
504
505#endif // THYRA_EPETRA_MODEL_EVALUATOR_HPP
506
507#if defined(Thyra_SHOW_DEPRECATED_WARNINGS)
508#ifdef __GNUC__
509#warning "The ThyraEpetraExtAdapters package is deprecated"
510#endif
511#endif
512
Concrete Adapter subclass that takes an EpetraExt::ModelEvaluator object and wraps it as a Thyra::Mod...
RCP< Teuchos::ParameterList > unsetParameterList()
RCP< Teuchos::ParameterList > getNonconstParameterList()
void setStateFunctionScalingVec(const RCP< const Epetra_Vector > &stateFunctionScalingVec)
Set the state function scaling vector s_f (see above).
RCP< const Epetra_Vector > getStateFunctionScalingVec() const
Get the state function scaling vector s_f (see above).
RCP< EpetraModelEvaluator > epetraModelEvaluator(const RCP< const EpetraExt::ModelEvaluator > &epetraModel, const RCP< LinearOpWithSolveFactoryBase< double > > &W_factory)
void initialize(const RCP< const EpetraExt::ModelEvaluator > &epetraModel, const RCP< LinearOpWithSolveFactoryBase< double > > &W_factory)
RCP< const VectorSpaceBase< double > > get_f_space() const
void setStateVariableScalingVec(const RCP< const Epetra_Vector > &stateVariableScalingVec)
Set the state variable scaling vector s_x (see above).
RCP< LinearOpBase< double > > create_W_op() const
ModelEvaluatorBase::InArgs< double > getUpperBounds() const
void setNominalValues(const ModelEvaluatorBase::InArgs< double > &nominalValues)
Set the nominal values.
RCP< PreconditionerBase< double > > create_W_prec() const
Returns null currently.
RCP< const LinearOpWithSolveFactoryBase< double > > get_W_factory() const
ModelEvaluatorBase::InArgs< double > getNominalValues() const
ModelEvaluatorBase::DerivativeSupport convert(const EpetraExt::ModelEvaluator::DerivativeSupport &derivativeSupport)
RCP< const Teuchos::ParameterList > getParameterList() const
void reportFinalPoint(const ModelEvaluatorBase::InArgs< double > &finalPoint, const bool wasSolved)
ModelEvaluatorBase::InArgs< double > createInArgs() const
void setParameterList(RCP< Teuchos::ParameterList > const &paramList)
RCP< const VectorSpaceBase< double > > get_g_space(int j) const
RCP< const EpetraExt::ModelEvaluator > getEpetraModel() const
EpetraExt::ModelEvaluator::EDerivativeMultiVectorOrientation convert(const ModelEvaluatorBase::EDerivativeMultiVectorOrientation &mvOrientation)
void uninitialize(RCP< const EpetraExt::ModelEvaluator > *epetraModel=NULL, RCP< LinearOpWithSolveFactoryBase< double > > *W_factory=NULL)
RCP< const VectorSpaceBase< double > > get_x_space() const
RCP< const Epetra_Vector > getStateVariableInvScalingVec() const
Get the state variable scaling vector s_x (see above).
ModelEvaluatorBase::InArgs< double > getLowerBounds() const
EpetraExt::ModelEvaluator::Derivative convert(const ModelEvaluatorBase::Derivative< double > &derivative, const RCP< const Epetra_Map > &fnc_map, const RCP< const Epetra_Map > &var_map)
RCP< const Epetra_Vector > getStateVariableScalingVec() const
Get the inverse state variable scaling vector inv_s_x (see above).
RCP< const Teuchos::ParameterList > getValidParameters() const
Teuchos::ArrayView< const std::string > get_g_names(int j) const
ModelEvaluatorBase::DerivativeProperties convert(const EpetraExt::ModelEvaluator::DerivativeProperties &derivativeProperties)
RCP< const Teuchos::Array< std::string > > get_p_names(int l) const
RCP< const VectorSpaceBase< double > > get_p_space(int l) const
const ModelEvaluatorBase::InArgs< double > & getFinalPoint() const
ModelEvaluatorBase::EDerivativeMultiVectorOrientation convert(const EpetraExt::ModelEvaluator::EDerivativeMultiVectorOrientation &mvOrientation)
Base class for all linear operators.
Factory interface for creating LinearOpWithSolveBase objects from compatible LinearOpBase objects.
Determines the forms of a general derivative that are supported.
Simple aggregate class that stores a derivative object as a general linear operator or as a multi-vec...
Concrete aggregate class for all input arguments computable by a ModelEvaluator subclass object.
Simple aggregate class that stores a derivative object as a general linear operator or as a multi-vec...
Concrete aggregate class for all output arguments computable by a ModelEvaluator subclass object.
Default base class for concrete model evaluators.
Simple public strict containing properties of a derivative object.