fix most pedantic compiler warnings in the basic infrastructure

i.e., using clang 3.8 to compile the test suite with the following
flags:

```
-Weverything
-Wno-documentation
-Wno-documentation-unknown-command
-Wno-c++98-compat
-Wno-c++98-compat-pedantic
-Wno-undef
-Wno-padded
-Wno-global-constructors
-Wno-exit-time-destructors
-Wno-weak-vtables
-Wno-float-equal
```

should not produce any warnings anymore. In my opinion the only flag
which would produce beneficial warnings is -Wdocumentation. This has
not been fixed in this patch because writing documentation is left for
another day (or, more likely, year).

note that this patch consists of a heavy dose of the OPM_UNUSED macro
and plenty of static_casts (to fix signedness issues). Fixing the
singedness issues were quite a nightmare and the fact that the Dune
API is quite inconsistent in that regard was not exactly helpful. :/

Finally this patch includes quite a few formatting changes (e.g., all
occurences of 'T &t' should be changed to `T& t`) and some fixes for
minor issues which I've found during the excercise.

I've made sure that all unit tests the test suite still pass
successfully and I've made sure that flow_ebos still works for Norne
and that it did not regress w.r.t. performance.

(Note that this patch does not fix compiler warnings triggered `ebos`
and `flow_ebos` but only those caused by the basic infrastructure or
the unit tests.)

v2: fix the warnings that occur if the dune-localfunctions module is
    not available. thanks to [at]atgeirr for testing.
v3: fix dune 2.3 build issue
This commit is contained in:
Andreas Lauser
2016-11-07 15:14:07 +01:00
parent 9b53930557
commit ec4b6c82dd
20 changed files with 610 additions and 435 deletions

View File

@@ -36,6 +36,7 @@
#include <opm/material/fluidmatrixinteractions/LinearMaterial.hpp>
#include <opm/material/fluidmatrixinteractions/EffToAbsLaw.hpp>
#include <opm/material/fluidmatrixinteractions/MaterialTraits.hpp>
#include <opm/material/common/Unused.hpp>
#include <dune/grid/yaspgrid.hh>
#include <dune/grid/io/file/dgfparser/dgfyasp.hh>
@@ -174,7 +175,7 @@ public:
/*!
* \copydoc Doxygen::defaultProblemConstructor
*/
RichardsLensProblem(Simulator &simulator)
RichardsLensProblem(Simulator& simulator)
: ParentType(simulator)
, pnRef_(1e5)
{
@@ -270,21 +271,21 @@ public:
* \copydoc FvBaseMultiPhaseProblem::temperature
*/
template <class Context>
Scalar temperature(const Context &context, unsigned spaceIdx, unsigned timeIdx) const
Scalar temperature(const Context& context, unsigned spaceIdx, unsigned timeIdx) const
{ return temperature(context.globalSpaceIndex(spaceIdx, timeIdx), timeIdx); }
Scalar temperature(unsigned globalSpaceIdx, unsigned timeIdx) const
Scalar temperature(unsigned OPM_UNUSED globalSpaceIdx, unsigned OPM_UNUSED timeIdx) const
{ return 273.15 + 10; } // -> 10°C
/*!
* \copydoc FvBaseMultiPhaseProblem::intrinsicPermeability
*/
template <class Context>
const DimMatrix &intrinsicPermeability(const Context &context,
const DimMatrix& intrinsicPermeability(const Context& context,
unsigned spaceIdx,
unsigned timeIdx) const
{
const GlobalPosition &pos = context.pos(spaceIdx, timeIdx);
const GlobalPosition& pos = context.pos(spaceIdx, timeIdx);
if (isInLens_(pos))
return lensK_;
return outerK_;
@@ -294,14 +295,16 @@ public:
* \copydoc FvBaseMultiPhaseProblem::porosity
*/
template <class Context>
Scalar porosity(const Context &context, unsigned spaceIdx, unsigned timeIdx) const
Scalar porosity(const Context& OPM_UNUSED context,
unsigned OPM_UNUSED spaceIdx,
unsigned OPM_UNUSED timeIdx) const
{ return 0.4; }
/*!
* \copydoc FvBaseMultiPhaseProblem::materialLawParams
*/
template <class Context>
const MaterialLawParams &materialLawParams(const Context &context,
const MaterialLawParams& materialLawParams(const Context& context,
unsigned spaceIdx,
unsigned timeIdx) const
{
@@ -309,7 +312,8 @@ public:
return materialLawParams(globalSpaceIdx, timeIdx);
}
const MaterialLawParams& materialLawParams(unsigned globalSpaceIdx, unsigned timeIdx) const
const MaterialLawParams& materialLawParams(unsigned globalSpaceIdx,
unsigned OPM_UNUSED timeIdx) const
{
if (dofIsInLens_[globalSpaceIdx])
return lensMaterialParams_;
@@ -322,14 +326,15 @@ public:
* \copydetails Doxygen::contextParams
*/
template <class Context>
Scalar referencePressure(const Context &context,
Scalar referencePressure(const Context& context,
unsigned spaceIdx,
unsigned timeIdx) const
{ return referencePressure(context.globalSpaceIndex(spaceIdx, timeIdx), timeIdx); }
// the Richards model does not have an element context available at all places
// where the reference pressure is required...
Scalar referencePressure(unsigned globalSpaceIdx, unsigned timeIdx) const
Scalar referencePressure(unsigned OPM_UNUSED globalSpaceIdx,
unsigned OPM_UNUSED timeIdx) const
{ return pnRef_; }
//! \}
@@ -343,15 +348,15 @@ public:
* \copydoc FvBaseProblem::boundary
*/
template <class Context>
void boundary(BoundaryRateVector &values,
const Context &context,
void boundary(BoundaryRateVector& values,
const Context& context,
unsigned spaceIdx,
unsigned timeIdx) const
{
const auto &pos = context.pos(spaceIdx, timeIdx);
const auto& pos = context.pos(spaceIdx, timeIdx);
if (onLeftBoundary_(pos) || onRightBoundary_(pos)) {
const auto &materialParams = this->materialLawParams(context, spaceIdx, timeIdx);
const auto& materialParams = this->materialLawParams(context, spaceIdx, timeIdx);
Scalar Sw = 0.0;
Opm::ImmiscibleFluidState<Scalar, FluidSystem> fs;
@@ -388,12 +393,12 @@ public:
* \copydoc FvBaseProblem::initial
*/
template <class Context>
void initial(PrimaryVariables &values,
const Context &context,
void initial(PrimaryVariables& values,
const Context& context,
unsigned spaceIdx,
unsigned timeIdx) const
{
const auto &materialParams = this->materialLawParams(context, spaceIdx, timeIdx);
const auto& materialParams = this->materialLawParams(context, spaceIdx, timeIdx);
Scalar Sw = 0.0;
Opm::ImmiscibleFluidState<Scalar, FluidSystem> fs;
@@ -412,35 +417,35 @@ public:
* everywhere.
*/
template <class Context>
void source(RateVector &rate,
const Context &context,
unsigned spaceIdx,
unsigned timeIdx) const
void source(RateVector& rate,
const Context& OPM_UNUSED context,
unsigned OPM_UNUSED spaceIdx,
unsigned OPM_UNUSED timeIdx) const
{ rate = Scalar(0.0); }
//! \}
private:
bool onLeftBoundary_(const GlobalPosition &pos) const
bool onLeftBoundary_(const GlobalPosition& pos) const
{ return pos[0] < this->boundingBoxMin()[0] + eps_; }
bool onRightBoundary_(const GlobalPosition &pos) const
bool onRightBoundary_(const GlobalPosition& pos) const
{ return pos[0] > this->boundingBoxMax()[0] - eps_; }
bool onLowerBoundary_(const GlobalPosition &pos) const
bool onLowerBoundary_(const GlobalPosition& pos) const
{ return pos[1] < this->boundingBoxMin()[1] + eps_; }
bool onUpperBoundary_(const GlobalPosition &pos) const
bool onUpperBoundary_(const GlobalPosition& pos) const
{ return pos[1] > this->boundingBoxMax()[1] - eps_; }
bool onInlet_(const GlobalPosition &pos) const
bool onInlet_(const GlobalPosition& pos) const
{
Scalar width = this->boundingBoxMax()[0] - this->boundingBoxMin()[0];
Scalar lambda = (this->boundingBoxMax()[0] - pos[0]) / width;
return onUpperBoundary_(pos) && 0.5 < lambda && lambda < 2.0 / 3.0;
}
bool isInLens_(const GlobalPosition &pos) const
bool isInLens_(const GlobalPosition& pos) const
{
for (unsigned i = 0; i < dimWorld; ++i) {
if (pos[i] < lensLowerLeft_[i] || pos[i] > lensUpperRight_[i])