43 #include "EpetraExt_DiagonalTransientModel.hpp"
44 #include "Epetra_Comm.h"
45 #include "Epetra_Map.h"
46 #include "Epetra_CrsGraph.h"
47 #include "Epetra_CrsMatrix.h"
48 #include "Epetra_LocalMap.h"
49 #include "Teuchos_ParameterList.hpp"
50 #include "Teuchos_StandardParameterEntryValidators.hpp"
51 #include "Teuchos_Assert.hpp"
52 #include "Teuchos_as.hpp"
61 const std::string Implicit_name =
"Implicit";
62 const bool Implicit_default =
true;
64 const std::string Gamma_min_name =
"Gamma_min";
65 const double Gamma_min_default = -0.9;
67 const std::string Gamma_max_name =
"Gamma_max";
68 const double Gamma_max_default = -0.01;
70 const std::string Coeff_s_name =
"Coeff_s";
71 const std::string Coeff_s_default =
"{ 0.0 }";
73 const std::string Gamma_fit_name =
"Gamma_fit";
74 const std::string Gamma_fit_default =
"Linear";
76 const std::string NumElements_name =
"NumElements";
77 const int NumElements_default = 1;
79 const std::string x0_name =
"x0";
80 const double x0_default = 10.0;
82 const std::string ExactSolutionAsResponse_name =
"Exact Solution as Response";
83 const bool ExactSolutionAsResponse_default =
false;
87 double evalR(
const double& t,
const double& gamma,
const double& s )
89 return (exp(gamma*t)*sin(s*t));
94 double d_evalR_d_s(
const double& t,
const double& gamma,
const double& s )
96 return (exp(gamma*t)*cos(s*t)*t);
102 const double& x,
const double& t,
const double& gamma,
const double& s
105 return ( gamma*x + evalR(t,gamma,s) );
111 const double& x0,
const double& t,
const double& gamma,
const double& s
115 return ( x0 * exp(gamma*t) );
116 return ( exp(gamma*t) * (x0 + (1.0/s) * ( 1.0 - cos(s*t) ) ) );
126 const double& t,
const double& gamma,
const double& s
131 return ( -exp(gamma*t)/(s*s) * ( 1.0 - sin(s*t)*(s*t) - cos(s*t) ) );
135 class UnsetParameterVector {
137 ~UnsetParameterVector()
140 *vec_ = Teuchos::null;
142 UnsetParameterVector(
143 const RCP<RCP<const Epetra_Vector> > &vec
148 void setVector(
const RCP<RCP<const Epetra_Vector> > &vec )
153 RCP<RCP<const Epetra_Vector> > vec_;
160 namespace EpetraExt {
167 Teuchos::RCP<Epetra_Comm>
const& epetra_comm
169 : epetra_comm_(epetra_comm.assert_not_null()),
170 implicit_(Implicit_default),
171 numElements_(NumElements_default),
172 gamma_min_(Gamma_min_default),
173 gamma_max_(Gamma_max_default),
174 coeff_s_(Teuchos::fromStringToArray<double>(Coeff_s_default)),
175 gamma_fit_(GAMMA_FIT_LINEAR),
177 exactSolutionAsResponse_(ExactSolutionAsResponse_default),
184 Teuchos::RCP<const Epetra_Vector>
191 Teuchos::RCP<const Epetra_Vector>
193 const double t,
const Epetra_Vector *coeff_s_p
196 set_coeff_s_p(Teuchos::rcp(coeff_s_p,
false));
197 UnsetParameterVector unsetParameterVector(Teuchos::rcp(&coeff_s_p_,
false));
198 Teuchos::RCP<Epetra_Vector>
199 x_star_ptr = Teuchos::rcp(
new Epetra_Vector(*epetra_map_,
false));
200 Epetra_Vector& x_star = *x_star_ptr;
201 Epetra_Vector& gamma = *gamma_;
202 int myN = x_star.MyLength();
203 for (
int i=0 ; i<myN ; ++i ) {
204 x_star[i] = x_exact( x0_, t, gamma[i], coeff_s(i) );
210 Teuchos::RCP<const Epetra_MultiVector>
212 const double t,
const Epetra_Vector *coeff_s_p
215 set_coeff_s_p(Teuchos::rcp(coeff_s_p,
false));
216 UnsetParameterVector unsetParameterVector(Teuchos::rcp(&coeff_s_p_,
false));
217 Teuchos::RCP<Epetra_MultiVector>
218 dxds_star_ptr = Teuchos::rcp(
new Epetra_MultiVector(*epetra_map_,np_,
false));
219 Epetra_MultiVector& dxds_star = *dxds_star_ptr;
220 dxds_star.PutScalar(0.0);
221 Epetra_Vector& gamma = *gamma_;
222 int myN = dxds_star.MyLength();
223 for (
int i=0 ; i<myN ; ++i ) {
224 const int coeff_s_idx_i = this->coeff_s_idx(i);
225 (*dxds_star(coeff_s_idx_i))[i] = dxds_exact( t, gamma[i], coeff_s(i) );
230 return (dxds_star_ptr);
238 Teuchos::RCP<Teuchos::ParameterList>
const& paramList
241 using Teuchos::get;
using Teuchos::getIntegralValue;
242 using Teuchos::getArrayFromStringParameter;
243 TEUCHOS_TEST_FOR_EXCEPT( is_null(paramList) );
245 isIntialized_ =
false;
246 paramList_ = paramList;
247 implicit_ = get<bool>(*paramList_,Implicit_name);
248 numElements_ = get<int>(*paramList_,NumElements_name);
249 gamma_min_ = get<double>(*paramList_,Gamma_min_name);
250 gamma_max_ = get<double>(*paramList_,Gamma_max_name);
251 coeff_s_ = getArrayFromStringParameter<double>(*paramList_,Coeff_s_name);
252 gamma_fit_ = getIntegralValue<EGammaFit>(*paramList_,Gamma_fit_name);
253 x0_ = get<double>(*paramList_,x0_name);
254 exactSolutionAsResponse_ = get<bool>(*paramList_,ExactSolutionAsResponse_name);
259 Teuchos::RCP<Teuchos::ParameterList>
266 Teuchos::RCP<Teuchos::ParameterList>
269 Teuchos::RCP<Teuchos::ParameterList> _paramList = paramList_;
270 paramList_ = Teuchos::null;
275 Teuchos::RCP<const Teuchos::ParameterList>
282 Teuchos::RCP<const Teuchos::ParameterList>
285 using Teuchos::RCP;
using Teuchos::ParameterList;
286 using Teuchos::tuple;
287 using Teuchos::setIntParameter;
using Teuchos::setDoubleParameter;
288 using Teuchos::setStringToIntegralParameter;
289 static RCP<const ParameterList> validPL;
290 if (is_null(validPL)) {
291 RCP<ParameterList> pl = Teuchos::parameterList();
292 pl->set(Implicit_name,
true);
294 Gamma_min_name, Gamma_min_default,
"",
298 Gamma_max_name, Gamma_max_default,
"",
301 setStringToIntegralParameter<EGammaFit>(
302 Gamma_fit_name, Gamma_fit_default,
"",
303 tuple<std::string>(
"Linear",
"Random"),
304 tuple<EGammaFit>(GAMMA_FIT_LINEAR,GAMMA_FIT_RANDOM),
308 NumElements_name, NumElements_default,
"",
312 x0_name, x0_default,
"",
315 pl->set( Coeff_s_name, Coeff_s_default );
316 pl->set( ExactSolutionAsResponse_name, ExactSolutionAsResponse_default );
326 Teuchos::RCP<const Epetra_Map>
333 Teuchos::RCP<const Epetra_Map>
340 Teuchos::RCP<const Epetra_Map>
344 TEUCHOS_ASSERT_IN_RANGE_UPPER_EXCLUSIVE( l, 0, Np_ );
350 Teuchos::RCP<const Teuchos::Array<std::string> >
354 TEUCHOS_ASSERT_IN_RANGE_UPPER_EXCLUSIVE( l, 0, Np_ );
360 Teuchos::RCP<const Epetra_Map>
364 TEUCHOS_ASSERT_IN_RANGE_UPPER_EXCLUSIVE( j, 0, Ng_ );
370 Teuchos::RCP<const Epetra_Vector>
377 Teuchos::RCP<const Epetra_Vector>
384 Teuchos::RCP<const Epetra_Vector>
388 TEUCHOS_ASSERT( l == 0 );
394 Teuchos::RCP<Epetra_Operator>
398 return Teuchos::rcp(
new Epetra_CrsMatrix(::Copy,*W_graph_));
399 return Teuchos::null;
403 EpetraExt::ModelEvaluator::InArgs
408 inArgs.setSupports(IN_ARG_t,
true);
409 inArgs.setSupports(IN_ARG_x,
true);
411 inArgs.setSupports(IN_ARG_x_dot,
true);
412 inArgs.setSupports(IN_ARG_alpha,
true);
413 inArgs.setSupports(IN_ARG_beta,
true);
419 EpetraExt::ModelEvaluator::OutArgs
422 OutArgsSetup outArgs;
423 outArgs.set_Np_Ng(Np_,Ng_);
424 outArgs.setSupports(OUT_ARG_f,
true);
426 outArgs.setSupports(OUT_ARG_W,
true);
427 outArgs.set_W_properties(
428 DerivativeProperties(
429 DERIV_LINEARITY_NONCONST
435 outArgs.setSupports(OUT_ARG_DfDp,0,DERIV_MV_BY_COL);
436 outArgs.set_DfDp_properties(
437 0,DerivativeProperties(
438 DERIV_LINEARITY_NONCONST,
439 DERIV_RANK_DEFICIENT,
443 if (exactSolutionAsResponse_) {
444 outArgs.setSupports(OUT_ARG_DgDp,0,0,DERIV_MV_BY_COL);
445 outArgs.set_DgDp_properties(
446 0,0,DerivativeProperties(
447 DERIV_LINEARITY_NONCONST,
448 DERIV_RANK_DEFICIENT,
458 const InArgs& inArgs,
const OutArgs& outArgs
465 using Teuchos::dyn_cast;
467 const Epetra_Vector &x = *(inArgs.get_x());
468 const double t = inArgs.get_t();
469 if (Np_) set_coeff_s_p(inArgs.get_p(0));
470 UnsetParameterVector unsetParameterVector(rcp(&coeff_s_p_,
false));
475 Epetra_Operator *W_out = ( implicit_ ? outArgs.get_W().get() : 0 );
476 Epetra_Vector *f_out = outArgs.get_f().get();
477 Epetra_MultiVector *DfDp_out = 0;
478 if (Np_) DfDp_out = get_DfDp_mv(0,outArgs).get();
479 Epetra_Vector *g_out = 0;
480 Epetra_MultiVector *DgDp_out = 0;
481 if (exactSolutionAsResponse_) {
482 g_out = outArgs.get_g(0).get();
483 DgDp_out = get_DgDp_mv(0,0,outArgs,DERIV_MV_BY_COL).get();
486 const Epetra_Vector &gamma = *gamma_;
488 int localNumElements = x.MyLength();
491 Epetra_Vector &f = *f_out;
493 const Epetra_Vector *x_dot_in = inArgs.get_x_dot().get();
495 const Epetra_Vector &x_dot = *x_dot_in;
496 for (
int i=0 ; i<localNumElements ; ++i )
497 f[i] = x_dot[i] - f_func(x[i],t,gamma[i],coeff_s(i));
500 for (
int i=0 ; i<localNumElements ; ++i )
501 f[i] = - f_func(x[i],t,gamma[i],coeff_s(i));
505 for (
int i=0 ; i<localNumElements ; ++i )
506 f[i] = f_func(x[i],t,gamma[i],coeff_s(i));
512 const double alpha = inArgs.get_alpha();
513 const double beta = inArgs.get_beta();
514 Epetra_CrsMatrix &crsW = dyn_cast<Epetra_CrsMatrix>(*W_out);
518 = epetra_comm_->MyPID()*localNumElements + epetra_map_->IndexBase();
519 for(
int i = 0; i < localNumElements; ++i ) {
520 values[0] = alpha - beta*gamma[i];
521 indices[0] = i + offset;
522 crsW.ReplaceGlobalValues(
532 Epetra_MultiVector &DfDp = *DfDp_out;
535 for(
int i = 0; i < localNumElements; ++i ) {
536 DfDp[coeff_s_idx(i)][i]
537 = - d_evalR_d_s(t,gamma[i],coeff_s(i));
541 for(
int i = 0; i < localNumElements; ++i ) {
542 DfDp[coeff_s_idx(i)][i]
543 = + d_evalR_d_s(t,gamma[i],coeff_s(i));
566 void DiagonalTransientModel::initialize()
576 const int numProcs = epetra_comm_->NumProc();
577 const int procRank = epetra_comm_->MyPID();
578 epetra_map_ = rcp(
new Epetra_Map( numElements_ * numProcs, 0, *epetra_comm_ ) );
584 gamma_ = Teuchos::rcp(
new Epetra_Vector(*epetra_map_));
585 Epetra_Vector &gamma = *gamma_;
587 case GAMMA_FIT_LINEAR: {
588 const int N = gamma.GlobalLength();
589 const double slope = (gamma_max_ - gamma_min_)/(N-1);
590 const int MyLength = gamma.MyLength();
592 gamma[0] = gamma_min_;
595 for (
int i=0 ; i<MyLength ; ++i )
597 gamma[i] = slope*(procRank*MyLength+i)+gamma_min_;
602 case GAMMA_FIT_RANDOM: {
603 unsigned int seed = Teuchos::ScalarTraits<int>::random()+10*procRank;
608 const double slope = (gamma_min_ - gamma_max_)/2.0;
609 const double tmp = (gamma_max_ + gamma_min_)/2.0;
610 int MyLength = gamma.MyLength();
611 for (
int i=0 ; i<MyLength ; ++i)
613 gamma[i] = slope*gamma[i] + tmp;
618 TEUCHOS_TEST_FOR_EXCEPT(
"Should never get here!");
626 np_ = coeff_s_.size();
629 rcp(
new Epetra_LocalMap(np_,0,*epetra_comm_) )
633 Teuchos::RCP<Teuchos::Array<std::string> >
634 names_p_0 = Teuchos::rcp(
new Teuchos::Array<std::string>());
635 for (
int i = 0; i < np_; ++i )
636 names_p_0->push_back(
"coeff_s("+Teuchos::toString(i)+
")");
637 names_p_.push_back(names_p_0);
640 coeff_s_idx_.clear();
641 const int num_func_per_coeff_s = numElements_ / np_;
642 const int num_func_per_coeff_s_rem = numElements_ % np_;
643 for (
int coeff_s_idx_i = 0; coeff_s_idx_i < np_; ++coeff_s_idx_i ) {
644 const int num_func_per_coeff_s_idx_i
645 = num_func_per_coeff_s
646 + ( coeff_s_idx_i < num_func_per_coeff_s_rem ? 1 : 0 );
648 int coeff_s_idx_i_j = 0;
649 coeff_s_idx_i_j < num_func_per_coeff_s_idx_i;
653 coeff_s_idx_.push_back(coeff_s_idx_i);
657 TEUCHOS_TEST_FOR_EXCEPT(
658 ( as<int>(coeff_s_idx_.size()) != numElements_ ) &&
"Internal programming error!" );
665 if (exactSolutionAsResponse_) {
669 rcp(
new Epetra_LocalMap(1,0,*epetra_comm_) )
690 const int offset = procRank*numElements_ + epetra_map_->IndexBase();
692 for(
int i = 0; i < numElements_; ++i ) {
693 indices[0] = i + offset;
694 W_graph_->InsertGlobalIndices(
701 W_graph_->FillComplete();
710 x_init_ = Teuchos::rcp(
new Epetra_Vector(*epetra_map_,
false));
711 x_init_->PutScalar(x0_);
716 x_dot_init_ = Teuchos::rcp(
new Epetra_Vector(*epetra_map_,
false));
718 inArgs.set_x(x_init_);
721 outArgs.set_f(x_dot_init_);
723 x_dot_init_->Scale(-1.0);
729 rcp(
new Epetra_Vector( ::Copy, *map_p_[0], &coeff_s_[0] ) )
735 void DiagonalTransientModel::set_coeff_s_p(
736 const Teuchos::RCP<const Epetra_Vector> &coeff_s_p
739 if (!is_null(coeff_s_p))
740 coeff_s_p_ = coeff_s_p;
746 void DiagonalTransientModel::unset_coeff_s_p()
const
751 as<int>(
get_p_map(0)->NumGlobalElements()) == as<int>(coeff_s_.size()) );
753 coeff_s_p_ = Teuchos::rcp(
757 const_cast<double*>(&coeff_s_[0])
772 Teuchos::RCP<EpetraExt::DiagonalTransientModel>
773 EpetraExt::diagonalTransientModel(
774 Teuchos::RCP<Epetra_Comm>
const& epetra_comm,
775 Teuchos::RCP<Teuchos::ParameterList>
const& paramList
779 Teuchos::rcp(
new DiagonalTransientModel(epetra_comm));
780 if (!is_null(paramList))
781 diagonalTransientModel->setParameterList(paramList);
OutArgs createOutArgs() const
Teuchos::RCP< const Epetra_Vector > getExactSolution(const double t, const Epetra_Vector *coeff_s_p=0) const
Return the exact solution as a function of time.
Teuchos::RCP< const Epetra_Vector > get_gamma() const
Return the model vector gamma,.
Teuchos::RCP< const Epetra_Map > get_p_map(int l) const
Teuchos::RCP< const Teuchos::ParameterList > getParameterList() const
DiagonalTransientModel(Teuchos::RCP< Epetra_Comm > const &epetra_comm)
Teuchos::RCP< const Epetra_Vector > get_x_init() const
Teuchos::RCP< const Epetra_Map > get_f_map() const
Teuchos::RCP< const Teuchos::Array< std::string > > get_p_names(int l) const
Teuchos::RCP< Teuchos::ParameterList > getNonconstParameterList()
Teuchos::RCP< Teuchos::ParameterList > unsetParameterList()
Teuchos::RCP< DiagonalTransientModel > diagonalTransientModel(Teuchos::RCP< Epetra_Comm > const &epetra_comm, Teuchos::RCP< Teuchos::ParameterList > const ¶mList=Teuchos::null)
Nonmember constructor.
Teuchos::RCP< const Epetra_MultiVector > getExactSensSolution(const double t, const Epetra_Vector *coeff_s_p=0) const
Return the exact sensitivity of x as a function of time.
Teuchos::RCP< const Epetra_Map > get_x_map() const
InArgs createInArgs() const
Teuchos::RCP< Epetra_Operator > create_W() const
void setParameterList(Teuchos::RCP< Teuchos::ParameterList > const ¶mList)
Teuchos::RCP< const Epetra_Vector > get_x_dot_init() const
Teuchos::RCP< const Epetra_Vector > get_p_init(int l) const
void evalModel(const InArgs &inArgs, const OutArgs &outArgs) const
Teuchos::RCP< const Teuchos::ParameterList > getValidParameters() const
Teuchos::RCP< const Epetra_Map > get_g_map(int j) const