|
Rythmos - Transient Integration for Differential Equations
Version of the Day
|
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., 51 Franklin St, Fifth Floor, Boston, MA 02110-1301 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 #include "Teuchos_ParameterListAcceptorDefaultBase.hpp" 00034 #include "Teuchos_VerboseObjectParameterListHelpers.hpp" 00035 #include "Rythmos_InterpolationBufferHelpers.hpp" 00036 00037 00038 namespace Rythmos { 00039 00040 00044 template<class Scalar> 00045 class PointwiseInterpolationBufferAppender 00046 : virtual public InterpolationBufferAppenderBase<Scalar>, 00047 virtual public Teuchos::ParameterListAcceptorDefaultBase 00048 { 00049 public: 00050 00052 typedef typename Teuchos::ScalarTraits<Scalar>::magnitudeType ScalarMag; 00053 00057 void append( 00058 const InterpolationBufferBase<Scalar>& interpBuffSource, 00059 const TimeRange<Scalar>& range, 00060 const Ptr<InterpolationBufferBase<Scalar> > &interpBuffSink 00061 ); 00062 00065 00067 void describe( 00068 Teuchos::FancyOStream &out, 00069 const Teuchos::EVerbosityLevel verbLevel 00070 ) const; 00071 00073 00076 00078 void setParameterList(const RCP<ParameterList> ¶mList); 00079 00081 RCP<const ParameterList> getValidParameters() const; 00082 00084 00085 }; 00086 00087 00088 00093 template<class Scalar> 00094 RCP<PointwiseInterpolationBufferAppender<Scalar> > 00095 pointwiseInterpolationBufferAppender() 00096 { 00097 return Teuchos::rcp(new PointwiseInterpolationBufferAppender<Scalar>); 00098 } 00099 00100 00101 // 00102 // Implementations 00103 // 00104 00105 00106 template<class Scalar> 00107 void PointwiseInterpolationBufferAppender<Scalar>::append( 00108 const InterpolationBufferBase<Scalar>& interpBuffSource, 00109 const TimeRange<Scalar>& appendRange, 00110 const Ptr<InterpolationBufferBase<Scalar> > &interpBuffSink 00111 ) 00112 { 00113 TEUCHOS_ASSERT( !is_null(interpBuffSink) ); 00114 #ifdef HAVE_RYTHMOS_DEBUG 00115 this->assertAppendPreconditions(interpBuffSource,appendRange,*interpBuffSink); 00116 #endif // HAVE_RYTHMOS_DEBUG 00117 00118 RCP<Teuchos::FancyOStream> out = this->getOStream(); 00119 Teuchos::OSTab ostab(out,1,"PointwiseInterpolationBufferAppender::append"); 00120 if ( Teuchos::as<int>(this->getVerbLevel()) >= Teuchos::as<int>(Teuchos::VERB_HIGH) ) { 00121 *out << "Interpolation Buffer source range = [" << interpBuffSource.getTimeRange().lower() << "," << 00122 interpBuffSource.getTimeRange().upper() << "]" << std::endl; 00123 *out << "Append range = [" << appendRange.lower() << "," << appendRange.upper() << "]" << std::endl; 00124 *out << "Interpolation Buffer sink range = [" << interpBuffSink->getTimeRange().lower() << "," << 00125 interpBuffSink->getTimeRange().upper() << "]" << std::endl; 00126 } 00127 // Set up appendRange correctly to be either (] or [): 00128 RCP<const TimeRange<Scalar> > correctedAppendRange = Teuchos::rcp(&appendRange,false); 00129 if (compareTimeValues<Scalar>(interpBuffSink->getTimeRange().upper(),appendRange.lower()) == 0) { 00130 // adding to end of buffer 00131 correctedAppendRange = Teuchos::rcp(new TimeRange_oc<Scalar>(appendRange)); 00132 if ( Teuchos::as<int>(this->getVerbLevel()) >= Teuchos::as<int>(Teuchos::VERB_HIGH) ) { 00133 *out << "Corrected append range = (" << correctedAppendRange->lower() << "," << 00134 correctedAppendRange->upper() << "]" << std::endl; 00135 } 00136 } 00137 else if (compareTimeValues<Scalar>(interpBuffSink->getTimeRange().lower(),appendRange.upper()) == 0) { 00138 // adding to beginning of buffer 00139 correctedAppendRange = Teuchos::rcp(new TimeRange_co<Scalar>(appendRange)); 00140 if ( Teuchos::as<int>(this->getVerbLevel()) >= Teuchos::as<int>(Teuchos::VERB_HIGH) ) { 00141 *out << "Corrected append range = [" << correctedAppendRange->lower() << "," << 00142 correctedAppendRange->upper() << ")" << std::endl; 00143 } 00144 } 00145 00146 Array<Scalar> time_vec_in; 00147 interpBuffSource.getNodes(&time_vec_in); 00148 00149 Array<Scalar> time_vec; 00150 selectPointsInTimeRange(time_vec_in,*correctedAppendRange,Teuchos::outArg(time_vec)); 00151 if ( Teuchos::as<int>(this->getVerbLevel()) >= Teuchos::as<int>(Teuchos::VERB_HIGH) ) { 00152 *out << "Selected points for appending to sink buffer: " << time_vec << std::endl; 00153 } 00154 00155 Array<RCP<const Thyra::VectorBase<Scalar> > > x_vec; 00156 Array<RCP<const Thyra::VectorBase<Scalar> > > xdot_vec; 00157 Array<ScalarMag> accuracy_vec; 00158 interpBuffSource.getPoints(time_vec, &x_vec, &xdot_vec, &accuracy_vec); 00159 00160 if ( Teuchos::as<int>(this->getVerbLevel()) >= Teuchos::as<int>(Teuchos::VERB_HIGH) ) { 00161 *out << "Sink buffer range before addPoints = [" << interpBuffSink->getTimeRange().lower() << "," << 00162 interpBuffSink->getTimeRange().upper() << "]" << std::endl; 00163 } 00164 00165 interpBuffSink->addPoints(time_vec, x_vec, xdot_vec); 00166 00167 if ( Teuchos::as<int>(this->getVerbLevel()) >= Teuchos::as<int>(Teuchos::VERB_HIGH) ) { 00168 *out << "Sink buffer range after addPoints = [" << interpBuffSink->getTimeRange().lower() << "," << 00169 interpBuffSink->getTimeRange().upper() << "]" << std::endl; 00170 } 00171 00172 } 00173 00174 00175 template<class Scalar> 00176 void PointwiseInterpolationBufferAppender<Scalar>::describe( 00177 Teuchos::FancyOStream &out, 00178 const Teuchos::EVerbosityLevel verbLevel 00179 ) const 00180 { 00181 using Teuchos::as; 00182 if ( 00183 (as<int>(verbLevel) == as<int>(Teuchos::VERB_DEFAULT)) 00184 || (as<int>(verbLevel) >= as<int>(Teuchos::VERB_LOW)) 00185 ) 00186 { 00187 out << this->description() << std::endl; 00188 } 00189 } 00190 00191 00192 template<class Scalar> 00193 void PointwiseInterpolationBufferAppender<Scalar>::setParameterList( 00194 const RCP<ParameterList> ¶mList 00195 ) 00196 { 00197 TEUCHOS_TEST_FOR_EXCEPT( is_null(paramList) ); 00198 paramList->validateParameters(*this->getValidParameters()); 00199 Teuchos::readVerboseObjectSublist(&*paramList,this); 00200 setMyParamList(paramList); 00201 } 00202 00203 00204 template<class Scalar> 00205 RCP<const ParameterList> 00206 PointwiseInterpolationBufferAppender<Scalar>::getValidParameters() const 00207 { 00208 static RCP<Teuchos::ParameterList> validPL; 00209 if (is_null(validPL)) { 00210 RCP<Teuchos::ParameterList> pl = Teuchos::parameterList(); 00211 Teuchos::setupVerboseObjectSublist(&*pl); 00212 validPL = pl; 00213 } 00214 return (validPL); 00215 } 00216 00217 00218 } // namespace Rythmos 00219 00220 00221 #endif //RYTHMOS_POINTWISE_INTERPOLATION_BUFFER_APPENDER_HPP
1.7.6.1