Rythmos_PointwiseInterpolationBufferAppender.hpp

00001 //@HEADER
00002 // ***********************************************************************
00003 //
00004 //                     Rythmos Package
00005 //                 Copyright (2006) Sandia Corporation
00006 //
00007 // Under terms of Contract DE-AC04-94AL85000, there is a non-exclusive
00008 // license for use of this work by or on behalf of the U.S. Government.
00009 //
00010 // This library is free software; you can redistribute it and/or modify
00011 // it under the terms of the GNU Lesser General Public License as
00012 // published by the Free Software Foundation; either version 2.1 of the
00013 // License, or (at your option) any later version.
00014 //
00015 // This library is distributed in the hope that it will be useful, but
00016 // WITHOUT ANY WARRANTY; without even the implied warranty of
00017 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
00018 // Lesser General Public License for more details.
00019 //
00020 // You should have received a copy of the GNU Lesser General Public
00021 // License along with this library; if not, write to the Free Software
00022 // Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA 02111-1307
00023 // USA
00024 // Questions? Contact Todd S. Coffey (tscoffe@sandia.gov)
00025 //
00026 // ***********************************************************************
00027 //@HEADER
00028 
00029 #ifndef RYTHMOS_POINTWISE_INTERPOLATION_BUFFER_APPENDER_HPP
00030 #define RYTHMOS_POINTWISE_INTERPOLATION_BUFFER_APPENDER_HPP
00031 
00032 #include "Rythmos_InterpolationBufferAppenderBase.hpp"
00033 
00034 
00035 namespace Rythmos {
00036 
00037 
00041 template<class Scalar>
00042 class PointwiseInterpolationBufferAppender
00043   : virtual public InterpolationBufferAppenderBase<Scalar>
00044 {
00045 public:
00046 
00048   typedef typename Teuchos::ScalarTraits<Scalar>::magnitudeType ScalarMag;
00049 
00053   void append(
00054     const InterpolationBufferBase<Scalar>& interpBuffSource, 
00055     const TimeRange<Scalar>& range,
00056     const Ptr<InterpolationBufferBase<Scalar> > &interpBuffSink 
00057     );
00058 
00061 
00063   void describe(
00064     Teuchos::FancyOStream &out,
00065     const Teuchos::EVerbosityLevel verbLevel
00066     ) const;
00067 
00069 
00072 
00074   void setParameterList(RCP<ParameterList> const& paramList);
00076   RCP<ParameterList> getNonconstParameterList();
00078   RCP<ParameterList> unsetParameterList();
00079 
00081 
00082 private:
00083 
00084   RCP<ParameterList> parameterList_;
00085 
00086 };
00087 
00088 
00089 
00094 template<class Scalar>
00095 RCP<PointwiseInterpolationBufferAppender<Scalar> >
00096 pointwiseInterpolationBufferAppender()
00097 {
00098   return Teuchos::rcp(new PointwiseInterpolationBufferAppender<Scalar>);
00099 }
00100 
00101 
00102 //
00103 // Implementations
00104 //
00105 
00106 
00107 template<class Scalar>
00108 void PointwiseInterpolationBufferAppender<Scalar>::append(
00109   const InterpolationBufferBase<Scalar>& interpBuffSource, 
00110   const TimeRange<Scalar>& range,
00111   const Ptr<InterpolationBufferBase<Scalar> > &interpBuffSink 
00112   ) 
00113 {
00114 #ifdef TEUCHOS_DEBUG
00115   this->assertAppendPreconditions(interpBuffSource,range,*interpBuffSink);
00116 #endif // TEUCHOS_DEBUG
00117 
00118   Array<Scalar> time_vec_in;
00119   interpBuffSource.getNodes(&time_vec_in);
00120 
00121   Array<Scalar> time_vec;
00122   selectPointsInTimeRange(&time_vec,time_vec_in,range);
00123   // 2007/12/05: rabrtl: ToDo: Make the output argument last!
00124 
00125   Array<RCP<const Thyra::VectorBase<Scalar> > > x_vec;
00126   Array<RCP<const Thyra::VectorBase<Scalar> > > xdot_vec;
00127   Array<ScalarMag> accuracy_vec;
00128   interpBuffSource.getPoints(time_vec, &x_vec, &xdot_vec, &accuracy_vec);
00129 
00130   interpBuffSink->addPoints(time_vec, x_vec, xdot_vec);
00131 
00132 }
00133 
00134 
00135 template<class Scalar>
00136 void PointwiseInterpolationBufferAppender<Scalar>::describe(
00137   Teuchos::FancyOStream &out,
00138   const Teuchos::EVerbosityLevel verbLevel
00139   ) const
00140 {
00141   using Teuchos::as;
00142   if (
00143     (as<int>(verbLevel) == as<int>(Teuchos::VERB_DEFAULT))
00144     || (as<int>(verbLevel) >= as<int>(Teuchos::VERB_LOW))
00145     )
00146   {
00147     out << this->description() << std::endl;
00148   }
00149 }
00150 
00151 
00152 template<class Scalar>
00153 void PointwiseInterpolationBufferAppender<Scalar>::setParameterList(
00154   RCP<ParameterList> const& paramList
00155   )
00156 {
00157   TEST_FOR_EXCEPTION(
00158     paramList==Teuchos::null, std::logic_error,
00159     "Error, paramList == Teuchos::null!\n" );
00160   parameterList_ = paramList;
00161 }
00162 
00163 
00164 template<class Scalar>
00165 RCP<ParameterList>
00166 PointwiseInterpolationBufferAppender<Scalar>::getNonconstParameterList()
00167 {
00168   return(parameterList_);
00169 }
00170 
00171 
00172 template<class Scalar>
00173 RCP<ParameterList>
00174 PointwiseInterpolationBufferAppender<Scalar>::unsetParameterList()
00175 {
00176   RCP<ParameterList> temp_param_list = parameterList_;
00177   parameterList_ = Teuchos::null;
00178   return(temp_param_list);
00179 }
00180 
00181 
00182 } // namespace Rythmos
00183 
00184 
00185 #endif //RYTHMOS_POINTWISE_INTERPOLATION_BUFFER_APPENDER_HPP
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends

Generated on Tue Oct 20 10:24:09 2009 for Rythmos - Transient Integration for Differential Equations by  doxygen 1.6.1