Makes the time step control parallel.

The only stage where parallelism changes the adaptive time
stepping is when some inner products on the saturation and
pressure are computed.
This commit makes this part parallel by added an additonal boost::any
parameter to the time stepping and the controller. Per default this
is empty. In a parallel run it contains a ParallelIstlInformation object
encapsulating the information about the parallelisation. This then used
to compute the parallel inner product.
This commit is contained in:
Markus Blatt 2015-05-22 20:51:59 +02:00
parent b6d5a3cad5
commit 384ef84556
4 changed files with 61 additions and 16 deletions

View File

@ -37,7 +37,11 @@ namespace Opm {
{
public:
//! \brief contructor taking parameter object
AdaptiveTimeStepping( const parameter::ParameterGroup& param );
//! \param param The parameter object
//! \param pinfo The information about the data distribution
//! and communication for a parallel run.
AdaptiveTimeStepping( const parameter::ParameterGroup& param,
const boost::any& pinfo=boost::any() );
/** \brief step method that acts like the solver::step method
in a sub cycle of time steps

View File

@ -33,7 +33,8 @@ namespace Opm {
// AdaptiveTimeStepping
//---------------------
AdaptiveTimeStepping::AdaptiveTimeStepping( const parameter::ParameterGroup& param )
AdaptiveTimeStepping::AdaptiveTimeStepping( const parameter::ParameterGroup& param,
const boost::any& parallel_information )
: timeStepControl_()
, restart_factor_( param.getDefault("solver.restartfactor", double(0.1) ) )
, growth_factor_( param.getDefault("solver.growthfactor", double(1.25) ) )
@ -51,13 +52,13 @@ namespace Opm {
const double tol = param.getDefault("timestep.control.tol", double(1e-3) );
if( control == "pid" ) {
timeStepControl_ = TimeStepControlType( new PIDTimeStepControl( tol ) );
timeStepControl_ = TimeStepControlType( new PIDTimeStepControl( tol, parallel_information ) );
}
else if ( control == "pid+iteration" )
{
const int iterations = param.getDefault("timestep.control.targetiteration", defaultTargetIterations );
const double maxgrowth = param.getDefault("timestep.control.maxgrowth", double(3.0) );
timeStepControl_ = TimeStepControlType( new PIDAndIterationCountTimeStepControl( iterations, tol, maxgrowth ) );
timeStepControl_ = TimeStepControlType( new PIDAndIterationCountTimeStepControl( iterations, tol, maxgrowth, parallel_information ) );
}
else if ( control == "iterationcount" )
{

View File

@ -16,7 +16,7 @@
You should have received a copy of the GNU General Public License
along with OPM. If not, see <http://www.gnu.org/licenses/>.
*/
#include <config.h>
#include <cassert>
#include <cmath>
#include <iostream>
@ -80,12 +80,14 @@ namespace Opm
//
////////////////////////////////////////////////////////
PIDTimeStepControl::PIDTimeStepControl( const double tol, const bool verbose )
PIDTimeStepControl::PIDTimeStepControl( const double tol, const boost::any& pinfo,
const bool verbose )
: p0_()
, sat0_()
, tol_( tol )
, errors_( 3, tol_ )
, verbose_( verbose )
, parallel_information_(pinfo)
{}
void PIDTimeStepControl::initialize( const SimulatorState& state )
@ -114,11 +116,13 @@ namespace Opm
// compute || u^n - u^n+1 ||
const double stateOld = euclidianNormSquared( p0_.begin(), p0_.end() ) +
euclidianNormSquared( sat0_.begin(), sat0_.end() );
euclidianNormSquared( sat0_.begin(), sat0_.end(),
state.numPhases() );
// compute || u^n+1 ||
const double stateNew = euclidianNormSquared( state.pressure().begin(), state.pressure().end() ) +
euclidianNormSquared( state.saturation().begin(), state.saturation().end() );
euclidianNormSquared( state.saturation().begin(), state.saturation().end(),
state.numPhases() );
// shift errors
for( int i=0; i<2; ++i ) {
@ -164,8 +168,9 @@ namespace Opm
PIDAndIterationCountTimeStepControl( const int target_iterations,
const double tol,
const double maxgrowth,
const boost::any& pinfo,
const bool verbose)
: BaseType( tol, verbose )
: BaseType( tol, pinfo, verbose )
, target_iterations_( target_iterations )
, maxgrowth_( maxgrowth )
{}

View File

@ -21,7 +21,10 @@
#include <vector>
#include <boost/any.hpp>
#include <boost/range/iterator_range.hpp>
#include <opm/core/simulator/TimeStepControlInterface.hpp>
#include <opm/core/linalg/ParallelIstlInformation.hpp>
namespace Opm
{
@ -73,8 +76,12 @@ namespace Opm
/// \brief constructor
/// \param tol tolerance for the relative changes of the numerical solution to be accepted
/// in one time step (default is 1e-3)
/// \paramm pinfo The information about the parallel information. Needed to
/// compute parallel scalarproducts.
/// \param verbose if true get some output (default = false)
PIDTimeStepControl( const double tol = 1e-3, const bool verbose = false );
PIDTimeStepControl( const double tol = 1e-3,
const boost::any& pinfo = boost::any(),
const bool verbose = false );
/// \brief \copydoc TimeStepControlInterface::initialize
void initialize( const SimulatorState& state );
@ -83,15 +90,38 @@ namespace Opm
double computeTimeStepSize( const double dt, const int /* iterations */, const SimulatorState& state ) const;
protected:
// return inner product for given container, here std::vector
template <class Iterator>
double euclidianNormSquared( Iterator it, const Iterator end ) const
double euclidianNormSquared( Iterator it, const Iterator end,
int num_components=1 ) const
{
double product = 0.0 ;
for( ; it != end; ++it ) {
product += ( *it * *it );
#if HAVE_MPI
if ( parallel_information_.type() == typeid(ParallelISTLInformation) )
{
const ParallelISTLInformation& info =
boost::any_cast<const ParallelISTLInformation&>(parallel_information_);
std::size_t size_per_component = (end - it) / num_components;
assert((end - it) == num_components * size_per_component);
double component_product = 0.0;
for( std::size_t i = 0; i < num_components; ++i )
{
auto component_container =
boost::make_iterator_range(it + i * size_per_component,
it + (i + 1) * size_per_component);
info.computeReduction(component_container,
Opm::Reduction::makeInnerProductFunctor<double>(),
component_product);
}
return component_product;
}
else
#endif
{
double product = 0.0 ;
for( ; it != end; ++it ) {
product += ( *it * *it );
}
return product;
}
return product;
}
protected:
@ -102,6 +132,8 @@ namespace Opm
mutable std::vector< double > errors_;
const bool verbose_;
private:
const boost::any parallel_information_;
};
///////////////////////////////////////////////////////////////////////////////////////////////////////////////
@ -119,10 +151,13 @@ namespace Opm
/// \param tol tolerance for the relative changes of the numerical solution to be accepted
/// in one time step (default is 1e-3)
// \param maxgrowth max growth factor for new time step in relation of old time step (default = 3.0)
/// \paramm pinfo The information about the parallel information. Needed to
/// compute parallel scalarproducts.
/// \param verbose if true get some output (default = false)
PIDAndIterationCountTimeStepControl( const int target_iterations = 20,
const double tol = 1e-3,
const double maxgrowth = 3.0,
const boost::any& = boost::any(),
const bool verbose = false);
/// \brief \copydoc TimeStepControlInterface::computeTimeStepSize