29 #ifndef RYTHMOS_POINTWISE_INTERPOLATION_BUFFER_APPENDER_HPP
30 #define RYTHMOS_POINTWISE_INTERPOLATION_BUFFER_APPENDER_HPP
32 #include "Rythmos_InterpolationBufferAppenderBase.hpp"
33 #include "Teuchos_ParameterListAcceptorDefaultBase.hpp"
34 #include "Teuchos_VerboseObjectParameterListHelpers.hpp"
35 #include "Rythmos_InterpolationBufferHelpers.hpp"
44 template<
class Scalar>
47 virtual public Teuchos::ParameterListAcceptorDefaultBase
52 typedef typename Teuchos::ScalarTraits<Scalar>::magnitudeType
ScalarMag;
68 Teuchos::FancyOStream &out,
69 const Teuchos::EVerbosityLevel verbLevel
93 template<
class Scalar>
94 RCP<PointwiseInterpolationBufferAppender<Scalar> >
106 template<
class Scalar>
113 TEUCHOS_ASSERT( !is_null(interpBuffSink) );
114 #ifdef HAVE_RYTHMOS_DEBUG
115 this->assertAppendPreconditions(interpBuffSource,appendRange,*interpBuffSink);
116 #endif // HAVE_RYTHMOS_DEBUG
118 RCP<Teuchos::FancyOStream> out = this->getOStream();
119 Teuchos::OSTab ostab(out,1,
"PointwiseInterpolationBufferAppender::append");
120 if ( Teuchos::as<int>(this->getVerbLevel()) >= Teuchos::as<int>(Teuchos::VERB_HIGH) ) {
121 *out <<
"Interpolation Buffer source range = [" << interpBuffSource.
getTimeRange().lower() <<
"," <<
122 interpBuffSource.
getTimeRange().upper() <<
"]" << std::endl;
123 *out <<
"Append range = [" << appendRange.
lower() <<
"," << appendRange.
upper() <<
"]" << std::endl;
124 *out <<
"Interpolation Buffer sink range = [" << interpBuffSink->getTimeRange().lower() <<
"," <<
125 interpBuffSink->getTimeRange().upper() <<
"]" << std::endl;
128 RCP<const TimeRange<Scalar> > correctedAppendRange = Teuchos::rcp(&appendRange,
false);
129 if (compareTimeValues<Scalar>(interpBuffSink->getTimeRange().upper(),appendRange.
lower()) == 0) {
131 correctedAppendRange = Teuchos::rcp(
new TimeRange_oc<Scalar>(appendRange));
132 if ( Teuchos::as<int>(this->getVerbLevel()) >= Teuchos::as<int>(Teuchos::VERB_HIGH) ) {
133 *out <<
"Corrected append range = (" << correctedAppendRange->lower() <<
"," <<
134 correctedAppendRange->upper() <<
"]" << std::endl;
137 else if (compareTimeValues<Scalar>(interpBuffSink->getTimeRange().lower(),appendRange.
upper()) == 0) {
139 correctedAppendRange = Teuchos::rcp(
new TimeRange_co<Scalar>(appendRange));
140 if ( Teuchos::as<int>(this->getVerbLevel()) >= Teuchos::as<int>(Teuchos::VERB_HIGH) ) {
141 *out <<
"Corrected append range = [" << correctedAppendRange->lower() <<
"," <<
142 correctedAppendRange->upper() <<
")" << std::endl;
146 Array<Scalar> time_vec_in;
147 interpBuffSource.
getNodes(&time_vec_in);
149 Array<Scalar> time_vec;
150 selectPointsInTimeRange(time_vec_in,*correctedAppendRange,Teuchos::outArg(time_vec));
151 if ( Teuchos::as<int>(this->getVerbLevel()) >= Teuchos::as<int>(Teuchos::VERB_HIGH) ) {
152 *out <<
"Selected points for appending to sink buffer: " << time_vec << std::endl;
155 Array<RCP<const Thyra::VectorBase<Scalar> > > x_vec;
156 Array<RCP<const Thyra::VectorBase<Scalar> > > xdot_vec;
157 Array<ScalarMag> accuracy_vec;
158 interpBuffSource.
getPoints(time_vec, &x_vec, &xdot_vec, &accuracy_vec);
160 if ( Teuchos::as<int>(this->getVerbLevel()) >= Teuchos::as<int>(Teuchos::VERB_HIGH) ) {
161 *out <<
"Sink buffer range before addPoints = [" << interpBuffSink->getTimeRange().lower() <<
"," <<
162 interpBuffSink->getTimeRange().upper() <<
"]" << std::endl;
165 interpBuffSink->addPoints(time_vec, x_vec, xdot_vec);
167 if ( Teuchos::as<int>(this->getVerbLevel()) >= Teuchos::as<int>(Teuchos::VERB_HIGH) ) {
168 *out <<
"Sink buffer range after addPoints = [" << interpBuffSink->getTimeRange().lower() <<
"," <<
169 interpBuffSink->getTimeRange().upper() <<
"]" << std::endl;
175 template<
class Scalar>
177 Teuchos::FancyOStream &out,
178 const Teuchos::EVerbosityLevel verbLevel
183 (as<int>(verbLevel) == as<int>(Teuchos::VERB_DEFAULT))
184 || (as<int>(verbLevel) >= as<int>(Teuchos::VERB_LOW))
187 out << this->description() << std::endl;
192 template<
class Scalar>
194 const RCP<ParameterList> ¶mList
197 TEUCHOS_TEST_FOR_EXCEPT( is_null(paramList) );
198 paramList->validateParameters(*this->getValidParameters());
199 Teuchos::readVerboseObjectSublist(&*paramList,
this);
200 setMyParamList(paramList);
204 template<
class Scalar>
205 RCP<const ParameterList>
208 static RCP<Teuchos::ParameterList> validPL;
209 if (is_null(validPL)) {
210 RCP<Teuchos::ParameterList> pl = Teuchos::parameterList();
211 Teuchos::setupVerboseObjectSublist(&*pl);
221 #endif //RYTHMOS_POINTWISE_INTERPOLATION_BUFFER_APPENDER_HPP
Concrete InterplationBufferAppender subclass that just transfers notes without any regard for accurac...
virtual void getPoints(const Array< Scalar > &time_vec, Array< RCP< const Thyra::VectorBase< Scalar > > > *x_vec, Array< RCP< const Thyra::VectorBase< Scalar > > > *xdot_vec, Array< ScalarMag > *accuracy_vec) const =0
Get values from the buffer at different time points.
virtual void getNodes(Array< Scalar > *time_vec) const =0
Get interpolation nodes.
void append(const InterpolationBufferBase< Scalar > &interpBuffSource, const TimeRange< Scalar > &range, const Ptr< InterpolationBufferBase< Scalar > > &interpBuffSink)
Concrete implementation that simply copies the nodal points between the interpolation buffers...
void setParameterList(const RCP< ParameterList > ¶mList)
RCP< PointwiseInterpolationBufferAppender< Scalar > > pointwiseInterpolationBufferAppender()
Nonmember constructor function.
Base class for an interpolation buffer.
Base class for strategy objects that append data from one InterplationBufferBase object to another...
RCP< const ParameterList > getValidParameters() const
void describe(Teuchos::FancyOStream &out, const Teuchos::EVerbosityLevel verbLevel) const
virtual TimeRange< Scalar > getTimeRange() const =0
Return the range of time values where interpolation calls can be performed.
Teuchos::ScalarTraits< Scalar >::magnitudeType ScalarMag