30 #ifndef Rythmos_INTEGRATOR_BUILDER_DEF_H
31 #define Rythmos_INTEGRATOR_BUILDER_DEF_H
34 #include "Rythmos_IntegratorBuilder_decl.hpp"
35 #include "Rythmos_IntegrationControlStrategyAcceptingIntegratorBase.hpp"
36 #include "Rythmos_SolverAcceptingStepperBase.hpp"
39 #include "Teuchos_as.hpp"
42 #include "Rythmos_DefaultIntegrator.hpp"
43 #include "Rythmos_SimpleIntegrationControlStrategy.hpp"
44 #include "Rythmos_RampingIntegrationControlStrategy.hpp"
45 #include "Rythmos_FixedStepControlStrategy.hpp"
46 #include "Rythmos_SimpleStepControlStrategy.hpp"
47 #include "Rythmos_FirstOrderErrorStepControlStrategy.hpp"
48 #include "Rythmos_ImplicitBDFStepperStepControl.hpp"
49 #include "Rythmos_ImplicitBDFStepperRampingStepControl.hpp"
50 #include "Rythmos_InterpolationBuffer.hpp"
51 #include "Rythmos_PointwiseInterpolationBufferAppender.hpp"
53 #include "Rythmos_ImplicitBDFStepperErrWtVecCalc.hpp"
54 #include "Rythmos_LinearInterpolator.hpp"
55 #include "Rythmos_HermiteInterpolator.hpp"
56 #include "Rythmos_CubicSplineInterpolator.hpp"
59 #include "Rythmos_ForwardSensitivityStepper.hpp"
65 static std::string integratorSettings_name =
"Integrator Settings";
66 static std::string integratorSettings_docs =
67 "These parameters are used directly in setting up the Integrator.";
68 static std::string integratorSelection_name =
"Integrator Selection";
69 static std::string integratorSelection_docs =
70 "Select the Integrator to be used.";
71 static std::string integrationControlSelection_name =
72 "Integration Control Strategy Selection";
73 static std::string integrationControlSelection_docs =
74 "Note that some settings conflict between step control and integration "
75 "control. In general, the integration control decides which steps will "
76 "be fixed or variable, not the stepper. When the integration control "
77 "decides to take variable steps, the step control is then responsible "
78 "for choosing appropriate step-sizes.";
79 static std::string stepperSettings_name =
"Stepper Settings";
80 static std::string stepperSettings_docs =
81 "This parameter list sets various parameters for the Stepper.";
82 static std::string stepperSelection_name =
"Stepper Selection";
83 static std::string stepperSelection_docs =
84 "Selects the Stepper for the time integration. It should be that "
85 "some time integrators can be accessed through different Steppers, "
86 "e.g., Backward Euler can be obtained through the `Backward Euler', "
87 "a first-order `Implicit BDF', or a one-stage `Implicit RK' Stepper."
88 "Special note for `Implicit RK' Stepper: If a fully implicit RK Butcher "
89 "tableau is chosen, then the stepper will not be fully initialized "
90 "unless a W factory object is set on the IntegratorBuilder through "
92 static std::string ForwardEulerStepper_name =
"Forward Euler";
93 static std::string ForwardEulerStepper_docs =
94 "This is the basic Forward Euler method: x_n = x_{n-1} + dt*x_dot_{n-1}";
95 static std::string BackwardEulerStepper_name =
"Backward Euler";
96 static std::string BackwardEulerStepper_docs =
97 "This is the basic Backward Euler method: x_n = x_{n-1} + dt*x_dot_n";
98 static std::string ImplicitBDFStepper_name =
"Implicit BDF";
99 static std::string ImplicitBDFStepper_docs =
100 "This Stepper provides a re-implementation of the algorithms in "
101 "the LLNL Sundials code IDA. This is an implicit BDF integrator "
102 "for DAEs which uses variable step-sizes and variable-orders "
103 "first through fourth.";
104 static std::string rkButcherTableauSelection_name =
105 "Runge Kutta Butcher Tableau Selection";
106 static std::string rkButcherTableauSelection_docs =
107 "Only the Explicit RK Stepper and the Implicit RK Stepper accept an "
108 "RK Butcher Tableau.";
109 static std::string ExplicitRKStepper_name =
"Explicit RK";
110 static std::string ExplicitRKStepper_docs =
111 "This Stepper has many explicit time-integrators using Runge-Kutta "
112 "formulation and the Butcher Tableau specification. See `"
113 +rkButcherTableauSelection_name+
"' ParameterList for available options.";
114 static std::string ImplicitRKStepper_name =
"Implicit RK";
115 static std::string ImplicitRKStepper_docs =
116 "This Stepper has many implicit time-integrators using Runge-Kutta "
117 "formulation and the Butcher Tableau specification. See `"
118 +rkButcherTableauSelection_name+
"' ParameterList for available options.";
119 static std::string stepControlSettings_name =
"Step Control Settings";
120 static std::string stepControlSettings_docs =
121 "Not all step control strategies are compatible with each stepper. "
122 "If the strategy has the name of a stepper in its name, then it only "
123 "works with that stepper.";
124 static std::string stepControlSelection_name =
125 "Step Control Strategy Selection";
126 static std::string stepControlSelection_docs =
127 "Used to select the Control Strategy for the stepper.";
128 static std::string errWtVecSelection_name =
129 "Error Weight Vector Calculator Selection";
130 static std::string errWtVecSelection_docs =
131 "Not all ErrWtVec calculators are compatible with each step control "
132 "strategy. If the calculator has the name of a stepper or another "
133 "step control strategy in its name, then it only works with that step "
135 static std::string interpolationBufferSettings_name =
136 "Interpolation Buffer Settings";
137 static std::string interpolationBufferSettings_docs =
138 "This parameter list sets various parameters for the InterpolationBuffer.";
139 static std::string interpolationBufferSelection_name =
140 "Trailing Interpolation Buffer Selection";
141 static std::string interpolationBufferSelection_docs =
142 "Used to select the Interpolation Buffer.";
143 static std::string interpolationBufferAppenderSelection_name =
144 "Interpolation Buffer Appender Selection";
145 static std::string interpolationBufferAppenderSelection_docs =
146 "Used to select the Interpolation Buffer Appender.";
147 static std::string initialTime_name =
"Initial Time";
148 static int initialTime_default = 0;
149 static std::string initialTime_docs =
150 "The initial time to start integration.";
151 static std::string finalTimeRythmos_name =
"Final Time";
152 static int finalTimeRythmos_default = 1;
153 static std::string finalTimeRythmos_docs =
"The final time to end integration.";
154 static std::string landOnFinalTime_name =
"Land On Final Time";
155 static bool landOnFinalTime_default =
true;
156 static std::string landOnFinalTime_docs =
157 "Exactly land on the final time; do not step past final time and "
159 static std::string interpolatorSelection_name =
"Interpolator Selection";
160 static std::string interpolatorSelection_docs =
161 "Choose the interpolator to use.";
162 static std::string stepperInterpolatorSelection_docs =
163 "Note all Steppers accept an interpolator. Currently, only the "
164 "BackwardEuler stepper does.";
167 static std::string integratorBuilder_name =
"Rythmos::Integrator";
168 static std::string integratorBuilderType_name =
"Integrator Type";
169 static std::string integrationControlBuilder_name =
170 "Rythmos::IntegrationControlStrategy";
171 static std::string integrationControlBuilderType_name =
172 "Integration Control Strategy Type";
173 static std::string stepControlBuilder_name =
174 "Rythmos::StepControlStrategy";
175 static std::string stepControlBuilderType_name =
176 "Step Control Strategy Type";
177 static std::string interpolationBufferBuilder_name =
178 "Rythmos::InterpolationBuffer";
179 static std::string interpolationBufferBuilderType_name =
180 "Interpolation Buffer Type";
181 static std::string interpolationBufferAppenderBuilder_name =
182 "Rythmos::InterpolationBufferAppender";
183 static std::string interpolationBufferAppenderBuilderType_name =
184 "Interpolation Buffer Appender Type";
185 static std::string errWtVecCalcBuilder_name =
"Rythmos::ErrWtVecCalc";
186 static std::string errWtVecCalcBuilderType_name =
187 "Error Weight Vector Calculator Type";
188 static std::string interpolatorBuilder_name =
"Rythmos::Interpolator";
189 static std::string interpolatorBuilderType_name =
"Interpolator Type";
192 static std::string defaultIntegrator_name =
"Default Integrator";
193 static std::string defaultIntegrator_docs =
194 "This Integrator will accept an IntergationControlStrategy, and "
195 "can have an IntegrationObserver. The client can specify the "
196 "maximum number of time steps allowed. The Integrator will loop "
197 "over the Stepper until it reaches the requested time. For each "
198 "step, the step size will be determined through a couple "
199 "mechanisms/filters. If an Integration Control Strategy has "
200 "been specified, the step size and the step type (fixed or "
201 "variable) will be determined by it. Otherwise the step size "
202 "will be set to the maximum real value and the step type will "
203 "be variable. Next if the step size is beyond the final time "
204 "and the `"+landOnFinalTime_name+
"' is specified, the step size is "
205 "adjusted to advance the state to the final time. The Stepper "
206 "is passed this step size and type to advance the state. The "
207 "DefaultIntegrator determines the step size and type taken by "
208 "the Stepper, and if the step has failed. If the "
209 "IntegrationControlStrategy handles failures, it can suggest "
210 "another step size and retry with the Stepper. Otherwise, the "
211 "Integrator will fall through with a failure. With a successful "
212 "step of the Stepper, the Integrator repeats the above until it "
213 "reaches the requested time. Multiple requested times can be "
214 "passed to the Integrator.";
215 static std::string simpleIntegrationControl_name =
216 "Simple Integration Control Strategy";
217 static std::string simpleIntegrationControl_docs =
218 "This Integration Control Strategy is meant to be simple with "
219 "very few parameters controlling it. Basically the client can "
220 "select fixed step type (the Stepper can only take the requested "
221 "step size) or variable step type (the Stepper can adjust the step "
222 "size to meet accuracy, order, or other criteria). For fixed step "
223 "type, the client can specify the step size and number of steps. "
224 "For variable step type, the client can set the maximum step size "
226 static std::string rampingIntegrationControl_name =
227 "Ramping Integration Control Strategy";
228 static std::string rampingIntegrationControl_docs =
229 "This Integration Control Strategy is very similar to `"
230 +simpleIntegrationControl_name+
"' except for handling an initial "
231 "constant-sized steps followed by a ramped-fixed-sized steps, "
232 "and finally variable- or fixed-sized steps. The client needs to "
233 "additionally set the initial step size and the maximum number of "
234 "step failures allowed.";
235 static std::string rampErrIntegrationControl_name =
236 "Ramp and Error Integration Control Strategy";
237 static std::string fixedStepControl_name =
"Fixed Step Control Strategy";
238 static std::string fixedStepControl_docs =
239 "This Step Control Strategy can be used for Steppers setup for "
240 "variable step type (a stepper that can adjust its step size based "
241 "on accuracy, order or other criteria), but would like to make fixed "
242 "step sizes or used fixed step size as its default.\n";
243 static std::string simpleStepControl_name =
"Simple Step Control Strategy";
244 static std::string simpleStepControl_docs =
245 "This Step Control Strategy starts with the initial step size, "
246 "and simply increases or decreases the step size by the "
247 "appropriate factor which is based on the change in the "
248 "solution relative to the specified relative and absolute "
249 "tolerances (|dx| < r*|x| + a) and if solution status from the "
250 "solver passes. Additionally the step size is bounded by the "
251 "miminum and maximum step size, and the stepper will fail if "
252 "the step size fails more than the specified value.";
253 static std::string implicitBDFStepControl_name =
254 "Implicit BDF Stepper Step Control Strategy";
255 static std::string implicitBDFStepControl_docs =
256 "This Step Control Strategy is specifically for use with the `"
257 +ImplicitBDFStepper_name+
"' Stepper. The parameters in this list "
258 "and sublist are directly related to those available in SUNDAILS/IDA. "
259 "See Hindmarsh, `The PVODE and IDA Algorithms', 2000 for more details. ";
260 static std::string implicitBDFStepperErrWtVecCalc_name =
261 "Implicit BDF Stepper Error Weight Vector Calculator";
262 static std::string implicitBDFStepperErrWtVecCalc_docs =
263 "This Error Weight Vector Calculator is specifically for use with the `"
264 +ImplicitBDFStepper_name+
"' Stepper.";
265 static std::string firstOrderErrorStepControl_name =
266 "First Order Error Step Control Strategy";
267 static std::string firstOrderErrorStepControl_docs =
268 "This Step Control Strategy produces a step size based on a first-order "
269 "predictor (Forward Euler) and a first-order solution (Backward Euler) by "
270 "by using a weight norm of the difference between the predicted and "
271 "solution. See Gresho and Sani, `Incompressible Flow and the Finite "
272 "Element Method', Vol. 1, 1998, p. 268.";
273 static std::string implicitBDFRampingStepControl_name =
274 "Implicit BDF Stepper Ramping Step Control Strategy";
275 static std::string implicitBDFRampingStepControl_docs =
276 "This Step Control Strategy is specifically for use with the `"
277 +ImplicitBDFStepper_name+
"' Stepper, and has a two-phase approach: "
278 "constant step sizes and followed by variable step sizes. The step "
279 "size is adjusted based on the WRMS, see "
280 +implicitBDFStepperErrWtVecCalc_name;
281 static std::string defaultInterpolationBuffer_name =
"Interpolation Buffer";
282 static std::string defaultInterpolationBuffer_docs =
283 "Sets parameters for the Interpolation Buffer.";
284 static std::string pointwiseInterpolationBufferAppender_name =
285 "Pointwise Interpolation Buffer Appender";
286 static std::string pointwiseInterpolationBufferAppender_docs =
287 "Appender that just transfers nodes without any regard for accuracy or "
291 static std::string linearInterpolator_name =
"Linear Interpolator";
292 static std::string linearInterpolator_docs =
293 "This provides a simple linear interpolation between time nodes.";
294 static std::string hermiteInterpolator_name =
"Hermite Interpolator";
295 static std::string hermiteInterpolator_docs =
296 "This provides a piecewise cubic Hermite interpolation on each interval "
297 "where the data is the solution and its time derivatives at the end "
298 "points of the interval. It will match 3rd degree polynomials exactly "
299 "with both function values and derivatives.";
300 static std::string cubicSplineInterpolator_name =
"Cubic Spline Interpolator";
301 static std::string cubicSplineInterpolator_docs =
302 "This provides a cubic spline interpolation between time nodes.";
310 template<
class Scalar>
313 this->initializeDefaults_();
317 template<
class Scalar>
323 template<
class Scalar>
326 const std::string &integratorName
329 integratorBuilder_->setObjectFactory(integratorFactory, integratorName);
330 validPL_ = Teuchos::null;
334 template<
class Scalar>
337 const std::string &integrationControlName
340 integrationControlBuilder_->setObjectFactory(integrationControlFactory,
341 integrationControlName);
342 validPL_ = Teuchos::null;
346 template<
class Scalar>
348 const RCP<StepperBuilder<Scalar> > &stepperBuilder
351 TEUCHOS_TEST_FOR_EXCEPT(is_null(stepperBuilder));
352 stepperBuilder_ = stepperBuilder;
353 validPL_ = Teuchos::null;
357 template<
class Scalar>
360 return stepperBuilder_;
364 template<
class Scalar>
366 const RCP<RKButcherTableauBuilder<Scalar> > & rkbtBuilder
369 TEUCHOS_TEST_FOR_EXCEPT(is_null(rkbtBuilder));
370 rkbtBuilder_ = rkbtBuilder;
371 validPL_ = Teuchos::null;
375 template<
class Scalar>
378 stepControlStrategyFactory,
379 const std::string &stepControlName
382 stepControlBuilder_->setObjectFactory(stepControlStrategyFactory,
384 validPL_ = Teuchos::null;
388 template<
class Scalar>
391 const std::string &interpolationBufferName
394 interpolationBufferBuilder_->setObjectFactory(interpolationBufferFactory, interpolationBufferName);
395 validPL_ = Teuchos::null;
399 template<
class Scalar>
402 const std::string &interpolationBufferAppenderName
405 interpolationBufferAppenderBuilder_->setObjectFactory(
406 interpolationBufferAppenderFactory, interpolationBufferAppenderName);
407 validPL_ = Teuchos::null;
411 template<
class Scalar>
413 const RCP<
const Teuchos::AbstractFactory<ErrWtVecCalcBase<Scalar> > > &
415 const std::string &errWtVecCalcFactoryName
418 errWtVecCalcBuilder_->setObjectFactory(errWtVecCalcFactory,
419 errWtVecCalcFactoryName);
420 validPL_ = Teuchos::null;
424 template<
class Scalar>
428 const std::string &interpolatorFactoryName
431 interpolatorBuilder_->setObjectFactory(interpolatorFactory,
432 interpolatorFactoryName);
433 validPL_ = Teuchos::null;
437 template<
class Scalar>
439 const RCP<Thyra::LinearOpWithSolveFactoryBase<Scalar> > &wFactoryObject
442 TEUCHOS_ASSERT( !is_null(wFactoryObject) );
443 wFactoryObject_ = wFactoryObject;
447 template<
class Scalar>
449 RCP<Teuchos::ParameterList>
const& paramList
452 TEUCHOS_TEST_FOR_EXCEPT(is_null(paramList));
453 paramList->validateParameters(*this->getValidParameters());
454 paramList_ = paramList;
458 template<
class Scalar>
459 RCP<const Teuchos::ParameterList>
462 if (is_null(validPL_)) {
463 RCP<ParameterList> pl = Teuchos::parameterList();
466 ParameterList& integratorSettingsPL =
467 pl->sublist(integratorSettings_name,
false,integratorSettings_docs);
470 integratorSettingsPL.set(initialTime_name,
471 Teuchos::as<Scalar>(initialTime_default),
474 integratorSettingsPL.set(finalTimeRythmos_name,
475 Teuchos::as<Scalar>(finalTimeRythmos_default),
476 finalTimeRythmos_docs);
478 integratorSettingsPL.set(landOnFinalTime_name,landOnFinalTime_default,
479 landOnFinalTime_docs);
481 ParameterList& integratorSelectionPL =
482 integratorSettingsPL.sublist(integratorSelection_name,
false,
483 integratorSelection_docs)
484 .disableRecursiveValidation();
486 integratorSelectionPL.sublist(defaultIntegrator_name,
false,
487 defaultIntegrator_docs)
488 .disableRecursiveValidation();
489 integratorSelectionPL
490 .setParameters(*(integratorBuilder_->getValidParameters()));
494 ParameterList& integrationControlSelectionPL =
495 pl->sublist(integrationControlSelection_name,
false,
496 integrationControlSelection_docs)
497 .disableRecursiveValidation();
499 integrationControlSelectionPL.sublist(simpleIntegrationControl_name,
false,
500 simpleIntegrationControl_docs)
501 .disableRecursiveValidation();
503 integrationControlSelectionPL.sublist(rampingIntegrationControl_name,
505 rampingIntegrationControl_docs)
506 .disableRecursiveValidation();
507 integrationControlSelectionPL
508 .setParameters(*(integrationControlBuilder_->getValidParameters()));
511 ParameterList& stepperSettingsPL = pl->sublist(stepperSettings_name,
false,
512 stepperSettings_docs);
515 ParameterList& stepperSelectionPL =
516 stepperSettingsPL.sublist(stepperSelection_name,
false,
517 stepperSelection_docs)
518 .disableRecursiveValidation();
520 stepperSelectionPL.sublist(ForwardEulerStepper_name,
false,
521 ForwardEulerStepper_docs)
522 .disableRecursiveValidation();
524 stepperSelectionPL.sublist(BackwardEulerStepper_name,
false,
525 BackwardEulerStepper_docs)
526 .disableRecursiveValidation();
528 stepperSelectionPL.sublist(ImplicitBDFStepper_name,
false,
529 ImplicitBDFStepper_docs)
530 .disableRecursiveValidation();
532 stepperSelectionPL.sublist(ExplicitRKStepper_name,
false,
533 ExplicitRKStepper_docs)
534 .disableRecursiveValidation();
536 stepperSelectionPL.sublist(ImplicitRKStepper_name,
false,
537 ImplicitRKStepper_docs)
538 .disableRecursiveValidation();
540 .setParameters(*(stepperBuilder_->getValidParameters()));
542 ParameterList& stepControlSettingsPL =
543 stepperSettingsPL.sublist(stepControlSettings_name,
false,
544 stepControlSettings_docs);
547 ParameterList& stepControlSelectionPL =
548 stepControlSettingsPL.sublist(stepControlSelection_name,
false,
549 stepControlSelection_docs)
550 .disableRecursiveValidation();
552 stepControlSelectionPL.sublist(fixedStepControl_name,
false,
553 fixedStepControl_docs)
554 .disableRecursiveValidation();
556 stepControlSelectionPL.sublist(simpleStepControl_name,
false,
557 simpleStepControl_docs)
558 .disableRecursiveValidation();
560 stepControlSelectionPL.sublist(firstOrderErrorStepControl_name,
false,
561 firstOrderErrorStepControl_docs)
562 .disableRecursiveValidation();
564 stepControlSelectionPL.sublist(implicitBDFStepControl_name,
false,
565 implicitBDFStepControl_docs)
566 .disableRecursiveValidation();
568 stepControlSelectionPL.sublist(implicitBDFRampingStepControl_name,
570 implicitBDFRampingStepControl_docs)
571 .disableRecursiveValidation();
572 stepControlSelectionPL
573 .setParameters(*(stepControlBuilder_->getValidParameters()));
576 ParameterList& errWtVecSelectionPL =
577 stepControlSettingsPL.sublist(errWtVecSelection_name,
false,
578 errWtVecSelection_docs)
579 .disableRecursiveValidation();
581 errWtVecSelectionPL.sublist(implicitBDFStepperErrWtVecCalc_name,
583 implicitBDFStepperErrWtVecCalc_docs)
584 .disableRecursiveValidation();
586 .setParameters(*(errWtVecCalcBuilder_->getValidParameters()));
589 ParameterList& interpolatorSelectionPL =
590 stepperSettingsPL.sublist(interpolatorSelection_name,
false,
591 stepperInterpolatorSelection_docs)
592 .disableRecursiveValidation();
594 interpolatorSelectionPL.sublist(linearInterpolator_name,
false,
595 linearInterpolator_docs)
596 .disableRecursiveValidation();
598 interpolatorSelectionPL.sublist(hermiteInterpolator_name,
false,
599 hermiteInterpolator_docs)
600 .disableRecursiveValidation();
602 interpolatorSelectionPL.sublist(cubicSplineInterpolator_name,
false,
603 cubicSplineInterpolator_docs)
604 .disableRecursiveValidation();
605 interpolatorSelectionPL
606 .setParameters(*(interpolatorBuilder_->getValidParameters()));
609 ParameterList& rkbtSelectionPL =
610 stepperSettingsPL.sublist(rkButcherTableauSelection_name,
false,
611 rkButcherTableauSelection_docs)
612 .disableRecursiveValidation();
613 rkbtSelectionPL.setParameters(*(rkbtBuilder_->getValidParameters()));
618 ParameterList& interpolationBufferSettingsPL =
619 pl->sublist(interpolationBufferSettings_name,
false,
620 interpolationBufferSettings_docs);
623 ParameterList& interpolationBufferSelectionPL =
624 interpolationBufferSettingsPL.sublist(interpolationBufferSelection_name,
626 interpolationBufferSelection_docs)
627 .disableRecursiveValidation();
629 interpolationBufferSelectionPL
630 .sublist(defaultInterpolationBuffer_name,
false,
631 defaultInterpolationBuffer_docs)
632 .disableRecursiveValidation();
633 interpolationBufferSelectionPL
634 .setParameters(*(interpolationBufferBuilder_->getValidParameters()));
636 ParameterList& interpolationBufferAppenderSelectionPL =
637 interpolationBufferSettingsPL
638 .sublist(interpolationBufferAppenderSelection_name,
false,
639 interpolationBufferAppenderSelection_docs)
640 .disableRecursiveValidation();
642 interpolationBufferAppenderSelectionPL
643 .sublist(pointwiseInterpolationBufferAppender_name,
false,
644 pointwiseInterpolationBufferAppender_docs)
645 .disableRecursiveValidation();
646 interpolationBufferAppenderSelectionPL
647 .setParameters(*(interpolationBufferAppenderBuilder_->getValidParameters()));
649 ParameterList& interpolatorSelectionPL =
650 interpolationBufferSettingsPL.sublist(interpolatorSelection_name,
false,
651 interpolatorSelection_docs)
652 .disableRecursiveValidation();
654 interpolatorSelectionPL.sublist(linearInterpolator_name,
false,
655 linearInterpolator_docs)
656 .disableRecursiveValidation();
658 interpolatorSelectionPL.sublist(hermiteInterpolator_name,
false,
659 hermiteInterpolator_docs)
660 .disableRecursiveValidation();
662 interpolatorSelectionPL.sublist(cubicSplineInterpolator_name,
false,
663 cubicSplineInterpolator_docs)
664 .disableRecursiveValidation();
665 interpolatorSelectionPL
666 .setParameters(*(interpolatorBuilder_->getValidParameters()));
677 template<
class Scalar>
684 template<
class Scalar>
687 RCP<ParameterList> pl = paramList_;
688 paramList_ = Teuchos::null;
693 template<
class Scalar>
711 template<
class Scalar>
712 RCP<IntegratorBase<Scalar> >
714 const RCP<
const Thyra::ModelEvaluator<Scalar> >& model,
715 const Thyra::ModelEvaluatorBase::InArgs<Scalar>& initialCondition,
716 const RCP<Thyra::NonlinearSolverBase<Scalar> >& nlSolver
719 TEUCHOS_TEST_FOR_EXCEPTION( is_null(model), std::logic_error,
720 "Error! IntegratorBuilder::create(...) The model passed in is null!"
722 TEUCHOS_TEST_FOR_EXCEPTION( is_null(paramList_), std::logic_error,
723 "Error! IntegratorBuilder::create(...) Please set a parameter "
724 "list on this class before calling create."
726 RCP<ParameterList> integratorSettingsPL = sublist(paramList_,
727 integratorSettings_name);
730 RCP<ParameterList> integratorSelectionPL = sublist(integratorSettingsPL,
731 integratorSelection_name);
732 integratorBuilder_->setParameterList(integratorSelectionPL);
733 RCP<IntegratorBase<Scalar> > integrator = integratorBuilder_->create();
734 TEUCHOS_TEST_FOR_EXCEPTION( is_null(integrator), std::logic_error,
735 "Error! IntegratorBuilder::create(...) The integrator came back "
736 "null from the ObjectBuilder!"
740 RCP<IntegrationControlStrategyAcceptingIntegratorBase<Scalar> >
743 if (!is_null(icsaIntegrator)) {
744 RCP<ParameterList> integrationControlSelectionPL =
745 sublist(paramList_,integrationControlSelection_name);
746 integrationControlBuilder_->setParameterList(integrationControlSelectionPL);
747 RCP<IntegrationControlStrategyBase<Scalar> > integrationControl =
748 integrationControlBuilder_->create();
749 if (!is_null(integrationControl)) {
750 icsaIntegrator->setIntegrationControlStrategy(integrationControl);
753 RCP<ParameterList> interpolationBufferSettingsPL =
754 sublist(paramList_,interpolationBufferSettings_name);
757 RCP<TrailingInterpolationBufferAcceptingIntegratorBase<Scalar> >
760 if (!is_null(tibaIntegrator)) {
761 RCP<ParameterList> interpolationBufferSelectionPL =
762 sublist(interpolationBufferSettingsPL,interpolationBufferSelection_name);
763 interpolationBufferBuilder_->setParameterList(interpolationBufferSelectionPL);
764 RCP<InterpolationBufferBase<Scalar> > ib =
765 interpolationBufferBuilder_->create();
768 RCP<InterpolatorAcceptingObjectBase<Scalar> > iaobIB =
770 if (!is_null(iaobIB)) {
771 RCP<ParameterList> interpolatorSelectionPL =
772 sublist(interpolationBufferSettingsPL,interpolatorSelection_name);
773 interpolatorBuilder_->setParameterList(interpolatorSelectionPL);
774 RCP<InterpolatorBase<Scalar> > interpolator =
775 interpolatorBuilder_->create();
776 if (!is_null(interpolator)) {
777 iaobIB->setInterpolator(interpolator);
780 tibaIntegrator->setTrailingInterpolationBuffer(ib);
785 RCP<InterpolationBufferAppenderAcceptingIntegratorBase<Scalar> > ibaaIntegrator =
787 if (!is_null(ibaaIntegrator)) {
788 RCP<ParameterList> interpolationBufferAppenderSelectionPL =
789 sublist(interpolationBufferSettingsPL,
790 interpolationBufferAppenderSelection_name);
791 interpolationBufferAppenderBuilder_->setParameterList(interpolationBufferAppenderSelectionPL);
792 RCP<InterpolationBufferAppenderBase<Scalar> > interpolationBufferAppender =
793 interpolationBufferAppenderBuilder_->create();
794 if (!is_null(interpolationBufferAppender)) {
795 ibaaIntegrator->setInterpolationBufferAppender(interpolationBufferAppender);
798 RCP<ParameterList> stepperSettingsPL =
799 sublist(paramList_,stepperSettings_name);
802 RCP<ParameterList> stepperSelectionPL = sublist(stepperSettingsPL,
803 stepperSelection_name);
804 stepperBuilder_->setParameterList(stepperSelectionPL);
805 RCP<StepperBase<Scalar> > stepper = stepperBuilder_->create();
806 TEUCHOS_TEST_FOR_EXCEPTION( is_null(stepper), std::logic_error,
807 "Error! IntegratorBuilder::create(...) The stepper came back "
808 "null from the StepperBuilder!");
811 RCP<ParameterList> stepControlSettingsPL =
812 sublist(stepperSettingsPL,stepControlSettings_name);
813 RCP<StepControlStrategyAcceptingStepperBase<Scalar> > scsaStepper =
815 if (!is_null(scsaStepper)) {
816 RCP<ParameterList> stepControlSelectionPL =
817 sublist(stepControlSettingsPL,stepControlSelection_name);
818 stepControlBuilder_->setParameterList(stepControlSelectionPL);
819 RCP<StepControlStrategyBase<Scalar> > stepControl =
820 stepControlBuilder_->create();
821 if (!is_null(stepControl)) {
823 RCP<ErrWtVecCalcAcceptingStepControlStrategyBase<Scalar> >
826 if (!is_null(ewvcaStepControl)) {
827 RCP<ParameterList> errWtVecSelectionPL =
828 sublist(stepControlSettingsPL,errWtVecSelection_name);
829 errWtVecCalcBuilder_->setParameterList(errWtVecSelectionPL);
830 RCP<ErrWtVecCalcBase<Scalar> > errWtVecCalc =
831 errWtVecCalcBuilder_->create();
832 if (!is_null(errWtVecCalc)) {
833 ewvcaStepControl->setErrWtVecCalc(errWtVecCalc);
836 scsaStepper->setStepControlStrategy(stepControl);
841 RCP<InterpolatorAcceptingObjectBase<Scalar> > iaobStepper =
844 if (!is_null(iaobStepper)) {
845 RCP<ParameterList> interpolatorSelectionPL =
846 sublist(stepperSettingsPL,interpolatorSelection_name);
847 interpolatorBuilder_->setParameterList(interpolatorSelectionPL);
848 RCP<InterpolatorBase<Scalar> > interpolator =
849 interpolatorBuilder_->create();
850 if (!is_null(interpolator)) {
851 iaobStepper->setInterpolator(interpolator);
856 RCP<RKButcherTableauAcceptingStepperBase<Scalar> > rkbtaStepper =
858 if (!is_null(rkbtaStepper)) {
859 RCP<ParameterList> rkButcherTableauSelectionPL =
860 sublist(stepperSettingsPL,rkButcherTableauSelection_name);
861 rkbtBuilder_->setParameterList(rkButcherTableauSelectionPL);
862 RCP<RKButcherTableauBase<Scalar> > rkbt = rkbtBuilder_->create();
863 TEUCHOS_TEST_FOR_EXCEPTION( is_null(rkbt), std::logic_error,
864 "Error! IntegratorBuilder::create(...) The Stepper accepts a "
865 "RK Butcher Tableau, but none were specified!"
867 rkbtaStepper->setRKButcherTableau(rkbt);
871 RCP<ImplicitRKStepper<Scalar> > irkStepper =
873 if (!is_null(irkStepper)) {
874 if (!is_null(wFactoryObject_)) {
875 irkStepper->set_W_factory(wFactoryObject_);
881 stepper->setModel(model);
883 stepper->setInitialCondition(initialCondition);
885 RCP<SolverAcceptingStepperBase<Scalar> > saStepper =
888 if(!is_null(saStepper)) {
889 TEUCHOS_TEST_FOR_EXCEPTION( is_null(nlSolver), std::logic_error,
890 "Error! IntegratorBuilder::create(...) The nonlinear solver passed "
891 "in is null and the stepper is implicit!"
893 saStepper->setSolver(nlSolver);
895 Scalar finalTimeRythmos = integratorSettingsPL->get<Scalar>(
896 finalTimeRythmos_name, Teuchos::as<Scalar>(finalTimeRythmos_default));
897 bool landOnFinalTime = integratorSettingsPL->get<
bool>(
898 landOnFinalTime_name, landOnFinalTime_default);
899 integrator->setStepper(stepper,finalTimeRythmos,landOnFinalTime);
903 template<
class Scalar>
907 using Teuchos::abstractFactoryStd;
910 integratorBuilder_ = Teuchos::objectBuilder<IntegratorBase<Scalar> >();
911 integratorBuilder_->setObjectName(integratorBuilder_name);
912 integratorBuilder_->setObjectTypeName(integratorBuilderType_name);
913 integratorBuilder_->setObjectFactory(
915 defaultIntegrator_name);
918 integrationControlBuilder_ =
919 Teuchos::objectBuilder<IntegrationControlStrategyBase<Scalar> >();
920 integrationControlBuilder_->setObjectName(integrationControlBuilder_name);
921 integrationControlBuilder_->setObjectTypeName(integrationControlBuilderType_name);
922 integrationControlBuilder_->setObjectFactory(
925 simpleIntegrationControl_name);
926 integrationControlBuilder_->setObjectFactory(
929 rampingIntegrationControl_name);
930 integrationControlBuilder_->setDefaultObject(
"None");
933 stepperBuilder_ = stepperBuilder<Scalar>();
936 rkbtBuilder_ = rKButcherTableauBuilder<Scalar>();
939 stepControlBuilder_ =
940 Teuchos::objectBuilder<StepControlStrategyBase<Scalar> >();
941 stepControlBuilder_->setObjectName(stepControlBuilder_name);
942 stepControlBuilder_->setObjectTypeName(stepControlBuilderType_name);
943 stepControlBuilder_->setObjectFactory(
945 FixedStepControlStrategy<Scalar> >(),
946 fixedStepControl_name);
947 stepControlBuilder_->setObjectFactory(
949 SimpleStepControlStrategy<Scalar> >(),
950 simpleStepControl_name);
951 stepControlBuilder_->setObjectFactory(
954 firstOrderErrorStepControl_name);
955 stepControlBuilder_->setObjectFactory(
957 ImplicitBDFStepperStepControl<Scalar> >(),
958 implicitBDFStepControl_name);
959 stepControlBuilder_->setObjectFactory(
962 implicitBDFRampingStepControl_name);
963 stepControlBuilder_->setDefaultObject(
"None");
966 interpolationBufferBuilder_ =
967 Teuchos::objectBuilder<InterpolationBufferBase<Scalar> >();
968 interpolationBufferBuilder_->setObjectName(interpolationBufferBuilder_name);
969 interpolationBufferBuilder_->setObjectTypeName(
970 interpolationBufferBuilderType_name);
971 interpolationBufferBuilder_->setObjectFactory(
974 defaultInterpolationBuffer_name);
975 interpolationBufferBuilder_->setDefaultObject(
"None");
978 interpolationBufferAppenderBuilder_ =
979 Teuchos::objectBuilder<InterpolationBufferAppenderBase<Scalar> >();
980 interpolationBufferAppenderBuilder_->setObjectName(
981 interpolationBufferAppenderBuilder_name);
982 interpolationBufferAppenderBuilder_->setObjectTypeName(
983 interpolationBufferAppenderBuilderType_name);
988 interpolationBufferAppenderBuilder_->setObjectFactory(
991 pointwiseInterpolationBufferAppender_name
993 interpolationBufferAppenderBuilder_->setDefaultObject(
"None");
996 errWtVecCalcBuilder_ = Teuchos::objectBuilder<ErrWtVecCalcBase<Scalar> >();
997 errWtVecCalcBuilder_->setObjectName(errWtVecCalcBuilder_name);
998 errWtVecCalcBuilder_->setObjectTypeName(errWtVecCalcBuilderType_name);
999 errWtVecCalcBuilder_->setObjectFactory(
1000 abstractFactoryStd< ErrWtVecCalcBase<Scalar>,
1001 ImplicitBDFStepperErrWtVecCalc<Scalar> >(),
1002 implicitBDFStepperErrWtVecCalc_name);
1003 errWtVecCalcBuilder_->setDefaultObject(
"None");
1006 interpolatorBuilder_ = Teuchos::objectBuilder<InterpolatorBase<Scalar> >();
1007 interpolatorBuilder_->setObjectName(interpolatorBuilder_name);
1008 interpolatorBuilder_->setObjectTypeName(interpolatorBuilderType_name);
1009 interpolatorBuilder_->setObjectFactory(
1012 linearInterpolator_name);
1013 interpolatorBuilder_->setObjectFactory(
1016 hermiteInterpolator_name);
1017 interpolatorBuilder_->setObjectFactory(
1020 cubicSplineInterpolator_name);
1021 interpolatorBuilder_->setDefaultObject(
"None");
1029 template<
class Scalar>
1030 Teuchos::RCP<Rythmos::IntegratorBuilder<Scalar> >
1031 Rythmos::integratorBuilder()
1033 return rcp(
new IntegratorBuilder<Scalar>);
1037 template<
class Scalar>
1038 Teuchos::RCP<Rythmos::IntegratorBuilder<Scalar> >
1039 Rythmos::integratorBuilder(
const RCP<ParameterList> ¶mList)
1041 const RCP<IntegratorBuilder<Scalar> > ib = integratorBuilder<Scalar>();
1042 ib->setParameterList(paramList);
1046 template<
class Scalar>
1047 Teuchos::RCP<Rythmos::IntegratorBase<Scalar> > Rythmos::createForwardSensitivityIntegrator(
1048 const RCP<
const Thyra::ModelEvaluator<Scalar> >& model,
1050 const Thyra::ModelEvaluatorBase::InArgs<Scalar>& model_ic,
1051 const RCP<Thyra::NonlinearSolverBase<Scalar> >& nlSolver,
1052 const RCP<Teuchos::ParameterList>& integratorBuilderPL
1055 RCP<IntegratorBuilder<Scalar> > ib = integratorBuilder<Scalar>(integratorBuilderPL);
1056 RCP<IntegratorBase<Scalar> > sensIntegrator = ib->create(model,model_ic,nlSolver);
1057 RCP<ForwardSensitivityStepper<Scalar> > stateAndSensStepper =
1058 forwardSensitivityStepper<Scalar>();
1059 stateAndSensStepper->initializeSyncedSteppers(
1060 model, p_index, model_ic, sensIntegrator->getNonconstStepper(), nlSolver
1062 typedef Thyra::ModelEvaluatorBase MEB;
1063 MEB::InArgs<Scalar> state_and_sens_ic =
1065 stateAndSensStepper->setInitialCondition(state_and_sens_ic);
1066 sensIntegrator->setStepper(stateAndSensStepper, sensIntegrator->getFwdTimeRange().upper());
1067 return sensIntegrator;
1078 #define RYTHMOS_INTEGRATOR_BUILDER_INSTANT(SCALAR) \
1080 template class IntegratorBuilder< SCALAR >; \
1082 template RCP<IntegratorBuilder< SCALAR > > \
1083 integratorBuilder(); \
1085 template RCP<IntegratorBuilder< SCALAR > > \
1086 integratorBuilder(const RCP<ParameterList> ¶List); \
1088 template RCP<IntegratorBase< SCALAR > > \
1089 createForwardSensitivityIntegrator( \
1090 const RCP<const Thyra::ModelEvaluator< SCALAR > >& model, \
1091 const int& p_index, \
1092 const Thyra::ModelEvaluatorBase::InArgs< SCALAR >& model_ic, \
1093 const RCP<Thyra::NonlinearSolverBase< SCALAR > >& nlSolver, \
1094 const RCP<ParameterList>& integratorBuilderPL \
1098 #endif //Rythmos_INTEGRATOR_BUILDER_DEF_H
void setRKButcherTableauBuilder(const RCP< RKButcherTableauBuilder< Scalar > > &rkbtBuilder)
Set the RK Butcher Tableau Builder object.
Concrete integrator builder class.
Mix-in interface for stepper objects that accept a step control strategy object to be used for evalua...
void setParameterList(const RCP< Teuchos::ParameterList > ¶mList)
void setIntegrationControlFactory(const RCP< const AbstractFactory< IntegrationControlStrategyBase< Scalar > > > &integrationControlFactory, const std::string &integrationControlName)
Set a new Integration Control Strategy factory object.
virtual ~IntegratorBuilder()
Concrete InterplationBufferAppender subclass that just transfers notes without any regard for accurac...
Base strategy class for interpolation functionality.
Abstract interface for time integrators.
RCP< const ParameterList > getParameterList() const
RCP< IntegratorBase< Scalar > > create(const RCP< const Thyra::ModelEvaluator< Scalar > > &model, const Thyra::ModelEvaluatorBase::InArgs< Scalar > &initialCondition, const RCP< Thyra::NonlinearSolverBase< Scalar > > &nlSolver) const
Create an fully formed integrator ready to go.
RCP< const Teuchos::ParameterList > getValidParameters() const
RCP< StepperBuilder< Scalar > > getStepperBuilder()
Get the Stepper Builder object.
Base class for strategy objects that control integration by selecting step sizes for a stepper...
Base class for strategy objects that control integration by selecting step sizes for a stepper...
Mix-in interface for integrator objects that accept an interpolationBufferAppender object to be used ...
Mix-in interface for integrator objects that accept an integration control strategy object to be used...
The member functions in the StepControlStrategyBase move you between these states in the following fa...
Mix-in interface for objects that accept an interpolator object.
Concrete implemenation of InterpolatorBase just just does simple linear interploation.
A concrete subclass for IntegratorBase that allows a good deal of customization.
Controls inital ramping at a fixed or incrementing time step size.
void setStepControlFactory(const RCP< const AbstractFactory< StepControlStrategyBase< Scalar > > > &stepControlStrategyFactory, const std::string &stepControlName)
Set a new Step Control Strategy factory object.
Thyra::ModelEvaluatorBase::InArgs< Scalar > createStateAndSensInitialCondition(const ForwardSensitivityStepper< Scalar > &fwdSensStepper, const Thyra::ModelEvaluatorBase::InArgs< Scalar > &state_ic, const RCP< const Thyra::MultiVectorBase< Scalar > > S_init=Teuchos::null, const RCP< const Thyra::MultiVectorBase< Scalar > > S_dot_init=Teuchos::null)
Set up default initial conditions for the state and sensitivity stepper with default zero initial con...
void setWFactoryObject(const RCP< Thyra::LinearOpWithSolveFactoryBase< Scalar > > &wFactoryObject)
Set a W factory object.
Base class for an interpolation buffer.
void setInterpolatorFactory(const RCP< const AbstractFactory< InterpolatorBase< Scalar > > > &interpolatorFactory, const std::string &interpolatorFactoryName)
Set an Interpolator factory object.
void setIntegratorFactory(const RCP< const AbstractFactory< IntegratorBase< Scalar > > > &integratorFactory, const std::string &integratorFactoryName)
Set a new Integrator factory object.
Concrete implemenation of InterpolatorBase that implements cubic spline interpolation.
concrete class for interpolation buffer functionality.
Base class for strategy objects that append data from one InterplationBufferBase object to another...
RCP< ParameterList > unsetParameterList()
Step Control Strategy for first-order time integration.
Mix-in interface for step control strategy objects that accept an external error weight calculation a...
Mix-in interface stepper objects that accept an RK Butcher Tableau.
RCP< ParameterList > getNonconstParameterList()
void setStepperBuilder(const RCP< StepperBuilder< Scalar > > &stepperBuilder)
Set the Stepper Builder object.
void setInterpolationBufferFactory(const RCP< const AbstractFactory< InterpolationBufferBase< Scalar > > > &interpolationBufferFactory, const std::string &interpolationBufferName)
Set an InterpolationBuffer factory object.
void setErrWtVecCalcFactory(const RCP< const AbstractFactory< ErrWtVecCalcBase< Scalar > > > &errWtVecCalcFactory, const std::string &errWtVecCalcFactoryName)
Set an ErrWtVecCalc factory object.
void setInterpolationBufferAppenderFactory(const RCP< const AbstractFactory< InterpolationBufferAppenderBase< Scalar > > > &interpolationBufferAppenderFactory, const std::string &interpolationBufferAppenderName)
Set an InterpolationBufferAppender factory object.
Mix-in interface all implicit stepper objects that accept a nonlinear solver to be used to compute th...
Mix-in interface for integrator objects that accept a trailing interpolation buffer object to be used...