Skip to content
Snippets Groups Projects
Commit 4835a9b2 authored by Katherina Baber's avatar Katherina Baber
Browse files

introduced property "Scaling" to scale model, set to 1 by default

git-svn-id: svn://svn.iws.uni-stuttgart.de/DUMUX/dumux/trunk@5836 2fb0f335-1f38-0410-981e-8018bf24f1b0
parent 83b1e9fc
No related branches found
No related tags found
No related merge requests found
...@@ -109,7 +109,12 @@ class OnePTwoCBoxModel : public BoxModel<TypeTag> ...@@ -109,7 +109,12 @@ class OnePTwoCBoxModel : public BoxModel<TypeTag>
static constexpr Scalar upwindAlpha = GET_PROP_VALUE(TypeTag, PTAG(UpwindAlpha)); static constexpr Scalar upwindAlpha = GET_PROP_VALUE(TypeTag, PTAG(UpwindAlpha));
typedef Dune::FieldMatrix<Scalar, dimWorld, dimWorld> Tensor; typedef Dune::FieldMatrix<Scalar, dimWorld, dimWorld> Tensor;
const Scalar scale_;
public: public:
OnePTwoCBoxModel():
scale_(GET_PROP_VALUE(TypeTag, PTAG(Scaling)))
{}
/*! /*!
* \brief \copybrief Dumux::BoxModel::addOutputVtkFields * \brief \copybrief Dumux::BoxModel::addOutputVtkFields
* *
...@@ -122,9 +127,8 @@ public: ...@@ -122,9 +127,8 @@ public:
void addOutputVtkFields(const SolutionVector &sol, void addOutputVtkFields(const SolutionVector &sol,
MultiWriter &writer) MultiWriter &writer)
{ {
Scalar scale = 1;//e3;
typedef Dune::BlockVector<Dune::FieldVector<Scalar, 1> > ScalarField; typedef Dune::BlockVector<Dune::FieldVector<Scalar, 1> > ScalarField;
std::cout<<"Scaling darcy "<< scale_<<std::endl;
// create the required scalar fields // create the required scalar fields
unsigned numVertices = this->problem_().gridView().size(dim); unsigned numVertices = this->problem_().gridView().size(dim);
ScalarField *pressure = writer.template createField<Scalar, 1>(numVertices); ScalarField *pressure = writer.template createField<Scalar, 1>(numVertices);
...@@ -187,14 +191,14 @@ public: ...@@ -187,14 +191,14 @@ public:
i, i,
false); false);
(*pressure)[globalIdx] = volVars.pressure()*scale; (*pressure)[globalIdx] = volVars.pressure()*scale_;
(*delp)[globalIdx] = volVars.pressure()*scale - 1e5; (*delp)[globalIdx] = volVars.pressure()*scale_ - 1e5;
(*moleFrac0)[globalIdx] = volVars.moleFrac(0); (*moleFrac0)[globalIdx] = volVars.moleFrac(0);
(*moleFrac1)[globalIdx] = volVars.moleFrac(1); (*moleFrac1)[globalIdx] = volVars.moleFrac(1);
(*massFrac0)[globalIdx] = volVars.massFrac(0); (*massFrac0)[globalIdx] = volVars.massFrac(0);
(*massFrac1)[globalIdx] = volVars.massFrac(1); (*massFrac1)[globalIdx] = volVars.massFrac(1);
(*rho)[globalIdx] = volVars.density()*scale*scale*scale; (*rho)[globalIdx] = volVars.density()*scale_*scale_*scale_;
(*mu)[globalIdx] = volVars.viscosity()*scale; (*mu)[globalIdx] = volVars.viscosity()*scale_;
(*delFrac)[globalIdx] = volVars.massFrac(1)-volVars.moleFrac(1); (*delFrac)[globalIdx] = volVars.massFrac(1)-volVars.moleFrac(1);
}; };
...@@ -295,16 +299,16 @@ public: ...@@ -295,16 +299,16 @@ public:
//use vertiacl faces for vx and horizontal faces for vy calculation //use vertiacl faces for vx and horizontal faces for vy calculation
(*velocityX)[i] /= boxSurface[i][0]; (*velocityX)[i] /= boxSurface[i][0];
(*velocityX)[i] /= scale; (*velocityX)[i] /= scale_;
if (dim >= 2) if (dim >= 2)
{ {
(*velocityY)[i] /= boxSurface[i][1]; (*velocityY)[i] /= boxSurface[i][1];
(*velocityY)[i] /= scale; (*velocityY)[i] /= scale_;
} }
if (dim == 3) if (dim == 3)
{ {
(*velocityZ)[i] /= boxSurface[i][2]; (*velocityZ)[i] /= boxSurface[i][2];
(*velocityZ)[i] /= scale; (*velocityZ)[i] /= scale_;
} }
} }
#endif #endif
......
...@@ -60,6 +60,7 @@ NEW_PROP_TAG(FluidSystem); //!< Type of the multi-component relations ...@@ -60,6 +60,7 @@ NEW_PROP_TAG(FluidSystem); //!< Type of the multi-component relations
NEW_PROP_TAG(UpwindAlpha); //!< The default value of the upwind parameter NEW_PROP_TAG(UpwindAlpha); //!< The default value of the upwind parameter
NEW_PROP_TAG(EnableGravity); //!< Returns whether gravity is considered in the problem NEW_PROP_TAG(EnableGravity); //!< Returns whether gravity is considered in the problem
NEW_PROP_TAG(UseMoles); //!Defines whether mole (true) or mass (false) fractions are used NEW_PROP_TAG(UseMoles); //!Defines whether mole (true) or mass (false) fractions are used
NEW_PROP_TAG(Scaling); //!Defines Scaling of the model
} }
// \} // \}
} }
......
...@@ -54,6 +54,7 @@ namespace Properties ...@@ -54,6 +54,7 @@ namespace Properties
SET_INT_PROP(BoxOnePTwoC, NumEq, 2); //!< set the number of equations to 2 SET_INT_PROP(BoxOnePTwoC, NumEq, 2); //!< set the number of equations to 2
SET_INT_PROP(BoxOnePTwoC, NumPhases, 1); //!< The number of phases in the 1p2c model is 1 SET_INT_PROP(BoxOnePTwoC, NumPhases, 1); //!< The number of phases in the 1p2c model is 1
SET_INT_PROP(BoxOnePTwoC, NumComponents, 2); //!< The number of components in the 1p2c model is 2 SET_INT_PROP(BoxOnePTwoC, NumComponents, 2); //!< The number of components in the 1p2c model is 2
SET_SCALAR_PROP(BoxOnePTwoC, Scaling, 1); //!< Scaling of the model is set to 1 by default
//! Use the 1p2c local residual function for the 1p2c model //! Use the 1p2c local residual function for the 1p2c model
SET_TYPE_PROP(BoxOnePTwoC, LocalResidual, OnePTwoCLocalResidual<TypeTag>); SET_TYPE_PROP(BoxOnePTwoC, LocalResidual, OnePTwoCLocalResidual<TypeTag>);
......
0% Loading or .
You are about to add 0 people to the discussion. Proceed with caution.
Finish editing this message first!
Please register or to comment