// -*- mode: C++; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4 -*-
// vi: set et ts=4 sw=4 sts=4:
/*****************************************************************************
 *   See the file COPYING for full copying permissions.                      *
 *                                                                           *
 *   This program is free software: you can redistribute it and/or modify    *
 *   it under the terms of the GNU General Public License as published by    *
 *   the Free Software Foundation, either version 2 of the License, or       *
 *   (at your option) any later version.                                     *
 *                                                                           *
 *   This program is distributed in the hope that it will be useful,         *
 *   but WITHOUT ANY WARRANTY; without even the implied warranty of          *
 *   MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the            *
 *   GNU General Public License for more details.                            *
 *                                                                           *
 *   You should have received a copy of the GNU General Public License       *
 *   along with this program.  If not, see <http://www.gnu.org/licenses/>.   *
 *****************************************************************************/
/*!
 * \file
 *
 * \brief Soil contamination problem where DNAPL infiltrates a fully
 *        water saturated medium.
 */


#ifndef DUMUX_NTDEX108_TUTORIAL_PROBLEM_COUPLED_HH    // guardian macro /*@\label{tutorial-coupled:guardian1}@*/
#define DUMUX_NTDEX108_TUTORIAL_PROBLEM_COUPLED_HH    // guardian macro /*@\label{tutorial-coupled:guardian2}@*/


#if HAVE_UG
#include <dune/grid/uggrid.hh>
#include <dune/grid/io/file/dgfparser/dgfug.hh>
#endif

#include <dune/grid/yaspgrid.hh>
//#include <dumux/io/simplexgridcreator.hh>
#include <dumux/io/cubegridcreator.hh>
#include <dune/grid/alugrid.hh>
 #include <dune/grid/io/file/dgfparser/dgfalu.hh>

#include "h2ohcube.hh"
#include "co2hcube.hh"
#include <dumux/implicit/2p/2pmodel.hh>
#include <dumux/implicit/common/implicitporousmediaproblem.hh>
#include <dumux/implicit/cellcentered/ccpropertydefaults.hh>

//#include <dumux/material/fluidsystems/h2on2fluidsystem.hh>

#include "ntdex108_tutorialspatialparams_coupled.hh"
#include <dumux/linear/amgbackend.hh>

namespace Dumux
{

template <class TypeTag>
class NtdEx108TutorialProblemCoupled;

namespace Properties
{
NEW_TYPE_TAG(NtdEx108TutorialProblemCoupled, INHERITS_FROM(CCTwoP, NtdEx108TutorialSpatialParamsCoupled));

SET_PROP(NtdEx108TutorialProblemCoupled, Problem) {typedef Dumux::NtdEx108TutorialProblemCoupled<TypeTag> type;};

// set the GridCreator property
//SET_TYPE_PROP(NtdEx108TutorialProblemCoupled, GridCreator, Dumux::CubeGridCreator<TypeTag>);

// Set the grid type
SET_TYPE_PROP(NtdEx108TutorialProblemCoupled, Grid, Dune::ALUGrid<3,3, Dune::cube, Dune::nonconforming>);
//SET_TYPE_PROP(NtdEx108TutorialProblemCoupled, Grid, Dune::ALUGrid<3,3, Dune::cube, Dune::conforming>);
//#if HAVE_UG
//SET_TYPE_PROP(NtdEx108TutorialProblemCoupled, Grid, Dune::UGGrid<2>);

//#else
//SET_TYPE_PROP(NtdEx108TutorialProblemCoupled, Grid, Dune::YaspGrid<3>);
//#endif

// Set the problem property
//SET_TYPE_PROP(NtdEx108TutorialProblemCoupled, Problem, Dumux::NtdEx108TutorialProblemCoupled<TypeTag>);

// Set the wetting phase
SET_PROP(NtdEx108TutorialProblemCoupled, WettingPhase)
{
private:
    typedef typename GET_PROP_TYPE(TypeTag, Scalar) Scalar;
public:
    typedef Dumux::LiquidPhase<Scalar, Dumux::H2ONTD<Scalar> > type;
};

// Set the non-wetting phase
SET_PROP(NtdEx108TutorialProblemCoupled, NonwettingPhase)
{
private:
    typedef typename GET_PROP_TYPE(TypeTag, Scalar) Scalar;
public:
    typedef Dumux::LiquidPhase<Scalar, Dumux::CO2NTD<Scalar> > type;
};

// Enable partial reassembly of the jacobian matrix?
SET_BOOL_PROP(NtdEx108TutorialProblemCoupled, ImplicitEnablePartialReassemble, true);

// Enable reuse of jacobian matrices?
SET_BOOL_PROP(NtdEx108TutorialProblemCoupled, ImplicitEnableJacobianRecycling, false);

// Write the solutions of individual newton iterations?
SET_BOOL_PROP(NtdEx108TutorialProblemCoupled, NewtonWriteConvergence, false);

// Use forward differences instead of central differences
SET_INT_PROP(NtdEx108TutorialProblemCoupled, ImplicitNumericDifferenceMethod, +1);

// Linear solver settings
//SET_TYPE_PROP(NtdEx108TutorialProblemCoupled, LinearSolver, Dumux::BoxBiCGStabILU0Solver<TypeTag> );
SET_TYPE_PROP(NtdEx108TutorialProblemCoupled, LinearSolver, Dumux::AMGBackend<TypeTag> );
SET_INT_PROP(NtdEx108TutorialProblemCoupled, LinearSolverVerbosity, 0);
SET_INT_PROP(NtdEx108TutorialProblemCoupled, LinearSolverPreconditionerIterations, 3);
SET_SCALAR_PROP(NtdEx108TutorialProblemCoupled, LinearSolverPreconditionerRelaxation, 1.0);

// Enable gravity
SET_BOOL_PROP(NtdEx108TutorialProblemCoupled, ProblemEnableGravity, true);
}

template <class TypeTag >
class NtdEx108TutorialProblemCoupled : public ImplicitPorousMediaProblem<TypeTag>
{
    typedef ImplicitPorousMediaProblem<TypeTag> ParentType;
    typedef typename GET_PROP_TYPE(TypeTag, GridView) GridView;
    typedef typename GET_PROP_TYPE(TypeTag, Indices) Indices;
    typedef typename GET_PROP_TYPE(TypeTag, FluidSystem) FluidSystem;
    typedef typename GET_PROP_TYPE(TypeTag, WettingPhase) WettingPhase;
    typedef typename GET_PROP_TYPE(TypeTag, NonwettingPhase) NonwettingPhase;

    enum {
        // primary variable indices
        pwIdx = Indices::pwIdx,
        snIdx = Indices::snIdx,

        // equation indices
        contiNEqIdx = Indices::contiNEqIdx,

        // phase indices
        wPhaseIdx = Indices::wPhaseIdx,
        nPhaseIdx = Indices::nPhaseIdx,

        // world dimension
        dim = GridView::dimension,
        dimWorld = GridView::dimensionworld
    };


    typedef typename GET_PROP_TYPE(TypeTag, PrimaryVariables) PrimaryVariables;
    typedef typename GET_PROP_TYPE(TypeTag, BoundaryTypes) BoundaryTypes;
    typedef typename GET_PROP_TYPE(TypeTag, TimeManager) TimeManager;
    typedef typename GET_PROP_TYPE(TypeTag, Scalar) Scalar;
    typedef typename GridView::template Codim<dim>::Entity Vertex;
    typedef typename GridView::Intersection Intersection;
    typedef Dune::FieldVector<Scalar, dimWorld> GlobalPosition;


    typedef typename GridView::template Codim<0>::Entity Element;
    typedef typename GridView::template Codim<0>::Iterator ElementIterator; // For permeability output_NTD
    typedef typename GET_PROP_TYPE(TypeTag, GridCreator) GridCreator;
    typedef typename GET_PROP_TYPE(TypeTag, Grid) Grid;
    typedef Dune::GridPtr<Grid> GridPointer;

    typedef typename GET_PROP_TYPE(TypeTag, ElementVolumeVariables) ElementVolumeVariables;
    typedef typename GET_PROP_TYPE(TypeTag, FVElementGeometry) FVElementGeometry;
    typedef typename GET_PROP_TYPE(TypeTag, VolumeVariables) VolumeVariables; // For permeability output_NTD
    enum { isBox = GET_PROP_VALUE(TypeTag, ImplicitIsBox) }; // For permeability output_NTD
    enum { dofCodim = isBox ? dim : 0 }; // For permeability output_NTD

public:
    /*!
     * \brief The constructor
     *
     * \param timeManager The time manager
     * \param gridView The grid view
     */
    NtdEx108TutorialProblemCoupled(TimeManager &timeManager,
                const GridView &gridView)
    : ParentType(timeManager, gridView)
    {
        eps_ = 3e-6;
        temperature_ = 273.15 + 10; // -> 20°C
        name_ = GET_RUNTIME_PARAM_FROM_GROUP(TypeTag, std::string, Problem, Name);
        Scalar outputTimeInterval_ = 1E6;
        outputTimeInterval_ = GET_RUNTIME_PARAM_FROM_GROUP(TypeTag, Scalar, Problem, OutputTimeInterval);
        this->timeManager().startNextEpisode(outputTimeInterval_);
        freqRestart_ = GET_RUNTIME_PARAM_FROM_GROUP(TypeTag, int, Output, FreqRestart);

        GridPointer *gridPtr = &GridCreator::gridPtr();
        Grid *grid_ = &GridCreator::grid();

        std::cout << "Process " << GridCreator::grid().leafGridView().comm().rank() + 1
            << " has " << GridCreator::grid().leafGridView().size(0)
            << " elements and " << GridCreator::grid().leafGridView().size(dim)
            << " nodes.\n";

        this->spatialParams().loaddata(gridPtr); // load data in process 0

        //gridPtr->loadBalance(); // grid load balancing
        
        GridCreator::gridPtr().loadBalance(); // NTD modification 19/03/2015
        
        
        this->spatialParams().loaddata(gridPtr); // distribute data in each process
     
        std::cout << "Process " << GridCreator::grid().leafGridView().comm().rank() + 1
                  << " has " << GridCreator::grid().leafGridView().size(0)
                  << " elements and " << GridCreator::grid().leafGridView().size(dim)
                  << " nodes.\n"; 
    }

    /*!
     * \name Problem parameters
     */
    // \{

    /*!
     * \brief Returns the problem name
     *
     * This is used as a prefix for files generated by the simulation.
     */
    const char *name() const
    {
        return name_.c_str();
    }

    /*!
     * \brief User defined output after the time integration
     *
     * Will be called diretly after the time integration.
     */
    void postTimeStep()
    {
        // Calculate storage terms
        PrimaryVariables storage;
        this->model().globalStorage(storage);

        // Write mass balance information for rank 0
        if (this->gridView().comm().rank() == 0) {
            std::cout<<"Storage: " << storage << std::endl;
        }
        if ( this->gridView().comm().rank() == 0){
            std::ofstream outfile;
            outfile.open("TimeStepAll.txt", std::ios_base::app);
            outfile << this->timeManager().time() + this->timeManager().timeStepSize()<<"\n";
            outfile.close();
        }

        if ( this->gridView().comm().rank() == 0){
            std::ofstream outfilestor;
            outfilestor.open("StorageAll.txt", std::ios_base::app);
            outfilestor << storage << "\n";
            outfilestor.close();
        }
    }

    /*!
     * \brief Returns the temperature \f$ K \f$
     *
     * This problem assumes a uniform temperature of 20 degrees Celsius.
     */
    Scalar temperature() const
    { return temperature_; }

    //! Returns true if the current solution should be written to disk
    //! as a VTK file

    bool shouldWriteOutput() const
    {
        return
            ((this->timeManager().timeStepIndex() > 0) &&
            (this->timeManager().episodeWillBeOver()) ) ||
            (this->timeManager().willBeFinished ()) ;
    }

    void episodeEnd()
    {
        if (!this->timeManager().willBeFinished())
        this->timeManager().startNextEpisode(GET_RUNTIME_PARAM_FROM_GROUP(TypeTag, Scalar, Problem, OutputTimeInterval));
    }

    /*!
     * \brief Returns the source term
     *
     * \param values Stores the source values for the conservation equations in
     *               \f$ [ \textnormal{unit of primary variable} / (m^\textrm{dim} \cdot s )] \f$
     * \param globalPos The global position
     */

    /* 
    void sourceAtPos(PrimaryVariables &values,
                const GlobalPosition &globalPos) const
    {
        Scalar DeltaX_= GET_RUNTIME_PARAM(TypeTag, Scalar, Grid.UpperRightX)/GET_RUNTIME_PARAM(TypeTag, Scalar, Grid.NumberOfCellsX);
        Scalar DeltaY_= GET_RUNTIME_PARAM(TypeTag, Scalar, Grid.UpperRightY)/GET_RUNTIME_PARAM(TypeTag, Scalar, Grid.NumberOfCellsY);
        Scalar DeltaZ_= GET_RUNTIME_PARAM(TypeTag, Scalar, Grid.UpperRightZ)/GET_RUNTIME_PARAM(TypeTag, Scalar, Grid.NumberOfCellsZ);
        Scalar GazInjectionRate_ = GET_RUNTIME_PARAM(TypeTag, Scalar, Problem.InjectionRate);
        //Scalar WaterProductionRate_ = GET_RUNTIME_PARAM(TypeTag, Scalar, Problem.ProductionRate);
        Scalar RightX_ = this->bBoxMax()[0];
        Scalar RightY_ = this->bBoxMax()[1];

        //if ( (globalPos[0] <= DeltaX_+eps_) && (globalPos[1] >= RightY_ - DeltaY_ - eps_)){
        if ( (globalPos[0] >= RightX_/2 - DeltaX_-eps_) && (globalPos[1] >= RightY_/2 - DeltaY_ - eps_) && (globalPos[0] <= RightX_/2 + DeltaX_+eps_) && (globalPos[1] <= RightY_/2 + DeltaY_ + eps_)){
            values[Indices::contiWEqIdx] = 0.0;
            values[Indices::contiNEqIdx]= GazInjectionRate_;}
        else{
            values[Indices::contiWEqIdx] = 0.0;
            values[Indices::contiNEqIdx]= 0.0;}
    }

    */
    void solDependentSource(PrimaryVariables &values,
                      const Element &element,
                      const FVElementGeometry &fvGeometry,
                      const int scvIdx,
                      const ElementVolumeVariables &elemVolVars) const
    {
        Scalar DeltaX_= GET_RUNTIME_PARAM(TypeTag, Scalar, Grid.DeltaX);
        Scalar DeltaY_= GET_RUNTIME_PARAM(TypeTag, Scalar, Grid.DeltaY);
        Scalar DeltaZ_= GET_RUNTIME_PARAM(TypeTag, Scalar, Grid.DeltaZ);
        Scalar RightX_ = this->bBoxMax()[0];
        Scalar RightY_ = this->bBoxMax()[1];
        Scalar WellX_ = GET_RUNTIME_PARAM(TypeTag, Scalar, SpatialParams.WellCoordinateX);
        Scalar WellY_ = GET_RUNTIME_PARAM(TypeTag, Scalar, SpatialParams.WellCoordinateY);
        Scalar WellR_ = GET_RUNTIME_PARAM(TypeTag, Scalar, SpatialParams.WellRadius);
        //Scalar MassDebit_ = GET_RUNTIME_PARAM(TypeTag, Scalar, Problem.MassDebit);
        Scalar AnisoK_ = GET_RUNTIME_PARAM(TypeTag, Scalar, SpatialParams.AnisoK);

        const double N_PI = (3.14159265358979323846);
        const GlobalPosition &globalPos = element.geometry().corner(scvIdx);

        Scalar PressureRef_ = GET_RUNTIME_PARAM(TypeTag, Scalar, Problem.PressureRef);
        typename GET_PROP_TYPE(TypeTag, FluidState) fluidState;
        fluidState.setTemperature(temperature_);
        fluidState.setPressure(FluidSystem::wPhaseIdx, PressureRef_);
        fluidState.setPressure(FluidSystem::nPhaseIdx, PressureRef_);

        Scalar densityW_ = FluidSystem::density(fluidState, FluidSystem::wPhaseIdx);
        Scalar densityN_ = FluidSystem::density(fluidState, FluidSystem::nPhaseIdx);
        Scalar viscosityW_ = FluidSystem::viscosity(fluidState, FluidSystem::wPhaseIdx);
        Scalar viscosityN_ = FluidSystem::viscosity(fluidState, FluidSystem::nPhaseIdx);


        Scalar PressureHole_ = PressureRef_ + (GET_RUNTIME_PARAM(TypeTag, Scalar, Grid.UpperRightZ)-globalPos[dim-1])*9.81*1050;

        //int myintsInj[] = {2,3,4,7,23,24,25,27,28,29,30};
        int myintsInj[] = {1,2,20,21,22,23,24,25,26,27,28,29,30,31,32,33,34,35,36,71,72,73,74,75,76,77,78,79,80,81,82,83};
        //int myintsInj[] = {1,2,3,21,22,23,24,25,26,27,28,29,30,31,32,33,34,35,36,37,71,72,73,74,75,76,77,78,79,80,81,82,83,84};

        std::vector<int> injectionCells_ (myintsInj, myintsInj + sizeof(myintsInj)/sizeof(int) );

        Scalar mobilityGas_ = 1.0/8E-5*710;
        Scalar Eff_WellR_ = 0.28 * pow(pow(AnisoK_,1.0/2.0) * pow(DeltaX_,2.0) + pow(1.0/AnisoK_,1.0/2.0) * pow(DeltaY_,2.0) , 1.0/2.0) / (pow(AnisoK_,1.0/4.0) +  pow(1.0/AnisoK_,1.0/4.0));

        //if ((globalPos[0] <= WellX_+eps_) && (globalPos[0] >= WellX_-DeltaX_-eps_) && (globalPos[1] <= WellY_ + eps_) && (globalPos[1] >= WellY_-DeltaY_ - eps_)
        //        && (std::find(injectionCells_.begin(), injectionCells_.end(),  round((globalPos[2] + DeltaZ_/2)/DeltaZ_)) != injectionCells_.end()) ) {
        //
        if ((globalPos[0] <= WellX_+eps_) && (globalPos[0] >= WellX_-DeltaX_/2-eps_) && (globalPos[1] <= WellY_ + eps_) && (globalPos[1] >= WellY_-DeltaY_/2 - eps_)
               && (std::find(injectionCells_.begin(), injectionCells_.end(),  round((globalPos[2] + DeltaZ_/2)/DeltaZ_)) != injectionCells_.end()) ) {

                    values[Indices::contiWEqIdx] = 0.0;
                    Scalar PermXCell_ = 5E-13;
                    Scalar WellIndex_=((2*N_PI*PermXCell_*DeltaZ_) / (log(Eff_WellR_ / WellR_)));
                    values[Indices::contiNEqIdx] = mobilityGas_*WellIndex_*(1.5*PressureHole_ - elemVolVars[scvIdx].pressure(nPhaseIdx))/fvGeometry.subContVol[scvIdx].volume;
                    //values[Indices::contiNEqIdx] = GET_RUNTIME_PARAM(TypeTag, Scalar, Problem.InjectionRate);
                    //std::cout<<mobilityGas_*WellIndex_*(1.5*PressureHole_-elemVolVars[scvIdx].pressure(nPhaseIdx))/fvGeometry.subContVol[scvIdx].volume<<"\n";
                    //std::cout<<" Eff_Well_Radius  = "<<Eff_WellR_<<"\n";
                    //std::cout<<" MobilityGas = "<<mobilityGas_<<"\n";
                    //std::cout<<" Volume = "<<fvGeometry.subContVol[scvIdx].volume<<"\n";
                    //std::cout<<" WellIndex_ = "<<WellIndex_<<"\n";
                    //std::cout<<" Pressure(nPhaseIdx) = "<<PressureHole_<<"\n";
                    //std::cout<<" Pressure(nPhaseIdx) = "<<elemVolVars[scvIdx].pressure(nPhaseIdx)<<"\n";
                    }
        else
            values = 0.0;
    }
    
    /*!
     * \name Boundary conditions
     */
    // \{

    /*!
     * \brief Specifies which kind of boundary condition should be
     *        used for which equation on a given boundary segment
     *
     * \param values Stores the value of the boundary type
     * \param globalPos The global position
     */
    /*
    void boundaryTypesAtPos(BoundaryTypes &values,
            const GlobalPosition &globalPos) const
    {
        if (onLeftBoundary_(globalPos) || onRightBoundary_(globalPos)) {
            values.setAllDirichlet();
        }
        else {
            values.setAllNeumann();
        }
    }
    */

   void boundaryTypesAtPos(BoundaryTypes &bcTypes, const GlobalPosition &globalPos) const
    {
        Scalar DeltaX_= GET_RUNTIME_PARAM(TypeTag, Scalar, Grid.UpperRightX)/GET_RUNTIME_PARAM(TypeTag, Scalar, Grid.NumberOfCellsX);
        Scalar DeltaY_= GET_RUNTIME_PARAM(TypeTag, Scalar, Grid.UpperRightY)/GET_RUNTIME_PARAM(TypeTag, Scalar, Grid.NumberOfCellsY);
        Scalar DeltaZ_= GET_RUNTIME_PARAM(TypeTag, Scalar, Grid.UpperRightZ)/GET_RUNTIME_PARAM(TypeTag, Scalar, Grid.NumberOfCellsZ);
        Scalar RightX_ = this->bBoxMax()[0];
        Scalar RightY_ = this->bBoxMax()[1];

        //if ( (globalPos[0] >= RightX_ - DeltaX_ - eps_) && (globalPos[1] <= DeltaY_ + eps_))
        Scalar topY  = this->bBoxMax()[1];
        if ( (globalPos[1] > topY - eps_) || (globalPos[1] < eps_))
            bcTypes.setAllDirichlet();
        else
            bcTypes.setAllNeumann();
    }


    /*!
     * \brief Evaluates the boundary conditions for a Dirichlet
     *        boundary segment
     *
     * \param values Stores the Dirichlet values for the conservation equations in
     *               \f$ [ \textnormal{unit of primary variable} ] \f$
     * \param globalPos The global position
     */
    /*
    void dirichletAtPos(PrimaryVariables &values,
                        const GlobalPosition &globalPos) const
    {
        typename GET_PROP_TYPE(TypeTag, FluidState) fluidState;
        fluidState.setTemperature(temperature_);
        fluidState.setPressure(FluidSystem::wPhaseIdx, 1e5);
        fluidState.setPressure(FluidSystem::nPhaseIdx, 1e5);

        Scalar densityW = FluidSystem::density(fluidState, FluidSystem::wPhaseIdx);

        Scalar height = this->bBoxMax()[1] - this->bBoxMin()[1];
        Scalar depth = this->bBoxMax()[1] - globalPos[1];
        Scalar alpha = 1 + 1.5/height;
        Scalar width = this->bBoxMax()[0] - this->bBoxMin()[0];
        Scalar factor = (width*alpha + (1.0 - alpha)*globalPos[0])/width;

        // hydrostatic pressure scaled by alpha
        values[pwIdx] = 1e5 - factor*densityW*this->gravity()[1]*depth;
        values[snIdx] = 0.0;
    }
    */

    void dirichletAtPos(PrimaryVariables &values, const GlobalPosition &globalPos) const

    {
        typename GET_PROP_TYPE(TypeTag, FluidState) fluidState;
        Scalar densityW = FluidSystem::density(fluidState, FluidSystem::wPhaseIdx);
        Scalar depth = this->bBoxMax()[dim-1] - globalPos[dim-1];
        values[Indices::pwIdx] = GET_RUNTIME_PARAM(TypeTag, Scalar, Problem.PressureRef) - densityW*this->gravity()[dim-1]*depth; // p=p0+drho*g*h hydrostatic pressure, p0 = 5 bar
        //values[Indices::pwIdx] = GET_RUNTIME_PARAM(TypeTag, Scalar, Problem.PressureRef);
        values[Indices::snIdx] = 0;
        //values[Indices::snIdx] = GET_RUNTIME_PARAM(TypeTag, Scalar, SpatialParams.ImposedSaturation); // oil saturation on bottom boundary
        //values[Indices::snIdx] = 0;
    }
    /*!
     * \brief Evaluate the boundary conditions for a neumann
     *        boundary segment.
     *
     * \param values Stores the Neumann values for the conservation equations in
     *               \f$ [ \textnormal{unit of conserved quantity} / (m^(dim-1) \cdot s )] \f$
     * \param globalPos The position of the integration point of the boundary segment.
     *
     * For this method, the \a values parameter stores the mass flux
     * in normal direction of each phase. Negative values mean influx.
     */
    void neumannAtPos(PrimaryVariables &values,
                      const GlobalPosition &globalPos) const
    {
        values = 0.0;
    }
    // \}

    /*!
     * \name Volume terms
     */
    // \{


    /*!
     * \brief Evaluates the initial values for a control volume
     *
     * \param values Stores the initial values for the conservation equations in
     *               \f$ [ \textnormal{unit of primary variables} ] \f$
     * \param globalPos The global position
     */
    void initialAtPos(PrimaryVariables &values,
                      const GlobalPosition &globalPos) const
    {
        typename GET_PROP_TYPE(TypeTag, FluidState) fluidState;
        fluidState.setTemperature(temperature_);
        fluidState.setPressure(FluidSystem::wPhaseIdx, GET_RUNTIME_PARAM(TypeTag, Scalar, Problem.PressureRef));
        fluidState.setPressure(FluidSystem::nPhaseIdx, GET_RUNTIME_PARAM(TypeTag, Scalar, Problem.PressureRef));

        Scalar densityW = FluidSystem::density(fluidState, FluidSystem::wPhaseIdx);

        Scalar depth = this->bBoxMax()[dim-1] - globalPos[dim-1];
    
        // hydrostatic pressure
        values[pwIdx] = GET_RUNTIME_PARAM(TypeTag, Scalar, Problem.PressureRef) - densityW*this->gravity()[dim-1]*depth;
        //values[pwIdx] = GET_RUNTIME_PARAM(TypeTag, Scalar, Problem.PressureRef);
        values[snIdx] = 0.0;
    }
    // \}

// Modification_NTD for displaying permeability output
    /*
     * \brief Append all quantities of interest which can be derived
     *        from the solution of the current time step to the VTK
     *        writer.
     */
     
    void addOutputVtkFields()
        {
        typedef Dune::BlockVector<Dune::FieldVector<double, 1> > ScalarField;
         

         // get the number of degrees of freedom
         unsigned numDofs = this->model().numDofs();
         unsigned numElements = this->gridView().size(0);

         //create required scalar fields
         ScalarField *Kxx = this->resultWriter().allocateManagedBuffer(numDofs);
         ScalarField *Kzz = this->resultWriter().allocateManagedBuffer(numDofs);
         ScalarField *boxVolume = this->resultWriter().allocateManagedBuffer(numDofs);
         (*boxVolume) = 0;

         //Fill the scalar fields with values
         ScalarField *rank = this->resultWriter().allocateManagedBuffer(numElements);

         FVElementGeometry fvGeometry;
         VolumeVariables volVars;

         ElementIterator eIt = this->gridView().template begin<0>();
         ElementIterator eEndIt = this->gridView().template end<0>();
         for (; eIt != eEndIt; ++eIt)
         {
             int eIdx = this->elementMapper().map(*eIt);
             (*rank)[eIdx] = this->gridView().comm().rank();
             fvGeometry.update(this->gridView(), *eIt);


             for (int scvIdx = 0; scvIdx < fvGeometry.numScv; ++scvIdx)
             {
                 int globalIdx = this->model().dofMapper().map(*eIt, scvIdx, dofCodim);
                 volVars.update(this->model().curSol()[globalIdx],
                                *this,
                                *eIt,
                                fvGeometry,
                                scvIdx,
                                false);
                 (*boxVolume)[globalIdx] += fvGeometry.subContVol[scvIdx].volume;
                 // (*Kxx)[globalIdx] = this->spatialParams().intrinsicPermeability(*eIt, fvGeometry, scvIdx)[1][1];
                 Dune::FieldMatrix<Scalar, dim, dim> perm = this->spatialParams().intrinsicPermeability(*eIt, fvGeometry, scvIdx);
                 (*Kxx)[globalIdx] = perm[1][1];
                 (*Kzz)[globalIdx] = perm[dim-1][dim-1];
             }
         }

         //pass the scalar fields to the vtkwriter
         this->resultWriter().attachDofData(*Kxx, "Kxx", isBox);
         this->resultWriter().attachDofData(*Kzz, "Kzz", isBox);
         this->resultWriter().attachDofData(*boxVolume, "boxVolume", isBox);
       } 
// End of modification NTD

private:
    unsigned freqRestart_;
    Scalar temperature_;
    Scalar eps_;
    std::string name_;
};
} //end namespace

#endif
