Rename AutoDiff::ForwardBlock -> Opm::AutoDiffBlock.

Also moved AutoDiffHelpers.hpp content to Opm namespace, and modified other
files as required by these two changes.
This commit is contained in:
Atgeirr Flø Rasmussen
2013-09-19 12:53:28 +02:00
parent e9b933bf4f
commit 85f79c0e84
15 changed files with 118 additions and 105 deletions

View File

@@ -26,14 +26,17 @@
#include <iostream>
namespace Opm
{
// -------------------- class HelperOps --------------------
/// Contains vectors and sparse matrices that represent subsets or
/// operations on (AD or regular) vectors of data.
struct HelperOps
{
typedef AutoDiff::ForwardBlock<double>::M M;
typedef AutoDiff::ForwardBlock<double>::V V;
typedef AutoDiffBlock<double>::M M;
typedef AutoDiffBlock<double>::V V;
/// A list of internal faces.
typedef Eigen::Array<int, Eigen::Dynamic, 1> IFaces;
@@ -202,7 +205,7 @@ namespace {
template <typename Scalar>
class UpwindSelector {
public:
typedef AutoDiff::ForwardBlock<Scalar> ADB;
typedef AutoDiffBlock<Scalar> ADB;
UpwindSelector(const UnstructuredGrid& g,
const HelperOps& h,
@@ -297,11 +300,11 @@ namespace {
/// Returns x(indices).
template <typename Scalar, class IntVec>
AutoDiff::ForwardBlock<Scalar>
subset(const AutoDiff::ForwardBlock<Scalar>& x,
AutoDiffBlock<Scalar>
subset(const AutoDiffBlock<Scalar>& x,
const IntVec& indices)
{
return ::constructSubsetSparseMatrix<Scalar>(x.value().size(), indices) * x;
return constructSubsetSparseMatrix<Scalar>(x.value().size(), indices) * x;
}
@@ -312,19 +315,19 @@ Eigen::Array<Scalar, Eigen::Dynamic, 1>
subset(const Eigen::Array<Scalar, Eigen::Dynamic, 1>& x,
const IntVec& indices)
{
return (::constructSubsetSparseMatrix<Scalar>(x.size(), indices) * x.matrix()).array();
return (constructSubsetSparseMatrix<Scalar>(x.size(), indices) * x.matrix()).array();
}
/// Returns v where v(indices) == x, v(!indices) == 0 and v.size() == n.
template <typename Scalar, class IntVec>
AutoDiff::ForwardBlock<Scalar>
superset(const AutoDiff::ForwardBlock<Scalar>& x,
AutoDiffBlock<Scalar>
superset(const AutoDiffBlock<Scalar>& x,
const IntVec& indices,
const int n)
{
return ::constructSupersetSparseMatrix<Scalar>(n, indices) * x;
return constructSupersetSparseMatrix<Scalar>(n, indices) * x;
}
@@ -336,7 +339,7 @@ superset(const Eigen::Array<Scalar, Eigen::Dynamic, 1>& x,
const IntVec& indices,
const int n)
{
return ::constructSupersetSparseMatrix<Scalar>(n, indices) * x.matrix();
return constructSupersetSparseMatrix<Scalar>(n, indices) * x.matrix();
}
@@ -345,10 +348,10 @@ superset(const Eigen::Array<Scalar, Eigen::Dynamic, 1>& x,
/// elements of d on the diagonal.
/// Need to mark this as inline since it is defined in a header and not a template.
inline
AutoDiff::ForwardBlock<double>::M
spdiag(const AutoDiff::ForwardBlock<double>::V& d)
AutoDiffBlock<double>::M
spdiag(const AutoDiffBlock<double>::V& d)
{
typedef AutoDiff::ForwardBlock<double>::M M;
typedef AutoDiffBlock<double>::M M;
const int n = d.size();
M mat(n, n);
@@ -367,7 +370,7 @@ spdiag(const AutoDiff::ForwardBlock<double>::V& d)
template <typename Scalar>
class Selector {
public:
typedef AutoDiff::ForwardBlock<Scalar> ADB;
typedef AutoDiffBlock<Scalar> ADB;
Selector(const typename ADB::V& selection_basis)
{
@@ -421,10 +424,10 @@ spdiag(const AutoDiff::ForwardBlock<double>::V& d)
/// Returns the input expression, but with all Jacobians collapsed to one.
inline
AutoDiff::ForwardBlock<double>
collapseJacs(const AutoDiff::ForwardBlock<double>& x)
AutoDiffBlock<double>
collapseJacs(const AutoDiffBlock<double>& x)
{
typedef AutoDiff::ForwardBlock<double> ADB;
typedef AutoDiffBlock<double> ADB;
const int nb = x.numBlocks();
typedef Eigen::Triplet<double> Tri;
int nnz = 0;
@@ -457,9 +460,9 @@ collapseJacs(const AutoDiff::ForwardBlock<double>& x)
/// Returns the vertical concatenation [ x; y ] of the inputs.
inline
AutoDiff::ForwardBlock<double>
vertcat(const AutoDiff::ForwardBlock<double>& x,
const AutoDiff::ForwardBlock<double>& y)
AutoDiffBlock<double>
vertcat(const AutoDiffBlock<double>& x,
const AutoDiffBlock<double>& y)
{
const int nx = x.size();
const int ny = y.size();
@@ -585,4 +588,7 @@ inline Eigen::ArrayXd sign (const Eigen::ArrayXd& x)
return retval;
}
} // namespace Opm
#endif // OPM_AUTODIFFHELPERS_HEADER_INCLUDED