Hi Romain,

in statements like "globalPos[idx]" or "this->bboxMax()[idx]," the "idx" refers to the coordinate direction, so in 2D "0" means the x- and "1" the y-direction. The x-axis is directed towards the right and the y-axis towards the top.

This means that in the "boundaryTypes" function, you have to replace
else if (globalPos[1] < eps_ ...
by
else if (globalPos[1] > this->bboxMax()[1] - eps_ ...

In the "neumann" function, you have to replace
Scalar right = this->bboxMax()[1];
if (globalPos[1] > right - eps_ ) {
by
Scalar right = this->bboxMax()[0];
if (globalPos[0] > right - eps_ ) {

Kind regards
Bernd

On 01/14/2013 03:11 PM, CHASSAGNE Romain wrote:
Dear Dumux users,
I am learning Dumux and then I set up a test case, but I am stuck with boundaries conditions, I failed to manage what I wanted.. I have attached an image of what I look for (MyGoal.jpg), I managed to do MyCurrentState.jpg so far. It means that the top and bottom boundaries conditions have to be the other way around (bottom to top and top to bottom). I guess there are something to do with the green lines, I tried different things but it let the boundaries conditions unchanged, which is odd.. and then I contact you. The piece (concerning boundaries conditions) of the code I have done so far to getMyCurrentState.jpg from "tutorialproblem_coupled.hh" :
void boundaryTypes(BoundaryTypes &bcTypes, const Vertex &vertex) const
    {
        const GlobalPosition &globalPos = vertex.geometry().center();
        // left boundary, injector1
        if (globalPos[0] < eps_ && globalPos[1] < 1000)
           bcTypes.setAllDirichlet();
// Well here is the problem, meant to be the top boundary condition, currently it is not. else if (globalPos[1] < eps_ && globalPos[0] > 2500 && globalPos[0]<3500)
          bcTypes.setAllDirichlet();
        else // neuman for the remaining boundaries
           bcTypes.setAllNeumann();
    }
    //! Evaluates the Dirichlet boundary conditions for a finite volume
    //! on the grid boundary. Here, the 'values' parameter stores
    //! primary variables.
    void dirichlet(PrimaryVariables &values, const Vertex &vertex) const
    {
        values[Indices::pwIdx] = 200.0e3; // 200 kPa = 2 bar
values[Indices::SnIdx] = 0.0; // 0 % oil saturation on left boundary
    }
    //! Evaluates the boundary conditions for a Neumann boundary
    //! segment. Here, the 'values' parameter stores the mass flux in
    //! [kg/(m^2 * s)] in normal direction of each phase. Negative
    //! values mean influx.
    void neumann(PrimaryVariables &values,
                 const Element &element,
                 const FVElementGeometry &fvGeometry,
                 const Intersection &is,
                 int scvIdx,
                 int boundaryFaceIdx) const
    {
        const GlobalPosition &globalPos =
fvGeometry.boundaryFace[boundaryFaceIdx].ipGlobal;
        Scalar right = this->bboxMax()[1];
// extraction of oil on the right boundary for approx. 1.e6 seconds
        if (globalPos[1] > right - eps_ ) {
            // oil outflux of 30 g/(m * s) on the right boundary.
            values[Indices::contiWEqIdx] = 0;
            values[Indices::contiNEqIdx] = 3e-2;
        } else {
            // no-flow on the remaining Neumann-boundaries.
            values[Indices::contiWEqIdx] = 0;
            values[Indices::contiNEqIdx] = 0;
        }
Thanks to the Dumux community
Cheers,
Romain.

__________________________

/Avant d'imprimer, pensez à l'environnement ! Please consider the environment before printing ! / /Ce message et toutes ses pièces jointes sont confidentiels et établis à l'intention exclusive de ses destinataires. Toute utilisation non conforme à sa destination, toute diffusion ou toute publication, totale ou partielle, est interdite, sauf autorisation expresse. IFP Energies nouvelles décline toute responsabilité au titre de ce message. This message and any attachments are confidential and intended solely for the addressees. Any unauthorised use or dissemination is prohibited. IFP Energies nouvelles should not be liable for this message./
__________________________



_______________________________________________
Dumux mailing list
[email protected]
https://listserv.uni-stuttgart.de/mailman/listinfo/dumux


--
_____________________________________________________________________

Bernd Flemisch                               phone: +49 711 685 69162
IWS, Universität Stuttgart                   fax:   +49 711 685 60430
Pfaffenwaldring 61                  email: [email protected]
D-70569 Stuttgart                  url: www.hydrosys.uni-stuttgart.de
_____________________________________________________________________

_______________________________________________
Dumux mailing list
[email protected]
https://listserv.uni-stuttgart.de/mailman/listinfo/dumux

Reply via email to