-- 
_________________________________________________

Timo Koch                                      phone: +49 711 685 64676
IWS, Universität Stuttgart                  fax:   +49 711 685 60430
Pfaffenwaldring 61         email: [email protected]
D-70569 Stuttgart             url: www.iws.uni-stuttgart.de/en/lh2/
_________________________________________________

> On 27. Mar 2020, at 14:36, Dmitry Pavlov <[email protected]> wrote:
> 
> Timo,
> 
> Thank you, my understanding has certainly improved.
> 
> To be sure: so I should implement neumann() to return a NumEqVector of two 
> numbers, the first being the flux of the wetting phase through the boundary, 
> and the second being the flux on the nonwetting phase though the same 
> boundary?
> 
> 
yes exactly!

Timo
> Best regards,
> 
> Dmitry
> On 26.03.2020 20:10, Timo Koch wrote:
>> 
>> -- 
>> _________________________________________________
>> 
>> Timo Koch                                      phone: +49 711 685 64676
>> IWS, Universität Stuttgart                  fax:   +49 711 685 60430
>> Pfaffenwaldring 61         email: [email protected] 
>> <mailto:[email protected]>
>> D-70569 Stuttgart             url: www.iws.uni-stuttgart.de/en/lh2/ 
>> <http://www.iws.uni-stuttgart.de/en/lh2/>
>> _________________________________________________
>> 
>>> On 26. Mar 2020, at 17:37, Dmitry Pavlov <[email protected] 
>>> <mailto:[email protected]>> wrote:
>>> 
>>> Dear Timo,
>>> 
>>> 
>>>> depends on what you understand under outflow.
>>> In the problem I am trying to model, I have one injection well and one 
>>> production well. The wells are essentially boundaries of a 2D rectangle. 
>>> Water is injected on one boundary, while oil+water mix is extracted on the 
>>> opposite one. I was thinking of setting on the production well a constant 
>>> pressure and zero relative saturation derivative (dSn/dx = 0, which of 
>>> course implies dSw/dx = 0). So this is a combination of Dirichlet for 
>>> pressure and Neumann for saturation. I am not quite sure that this is the 
>>> best way to model this sort of thing. If you have your thoughts, I would 
>>> appreciate them. I do not think that I can fix the saturation on the 
>>> production well, because I am studying a nonstationary situation in which 
>>> water gradually pushes out oil.
>>> 
>>> 
>>>> For Neumann what we expect from the function “neumann” in the problem 
>>>> class that it returns a boundary flux for both equations (NumEqVector). By 
>>>> “flux" I mean the entire term in the divergence.
>>> 
>>> Sorry, I am having trouble understanding how that can be applied to the 
>>> term with the saturation. But since I currently want to set it to zero, 
>>> probably it does not matter.
>> Hi Dmitry,
>> 
>> you have to return the normal flux term or what is left from it after 
>> inserting gradS*n = 0. That depends on which equation you are looking at, 
>> and what primary variables you have.
>> The flux term in the assembled equation (as you can see in the doc) contains 
>> grad(p_w) or grad(p_n). If you want to set something for S_w or S_n you 
>> first have to figure out (by using p_c = p_n - p_w, and S_w = S_w(p_c)) what 
>> the term looks like in terms of the saturation or its spatial derivatives.
>> Then insert your boundary conditions (pen&paper), and the left-over part has 
>> to be assembled as the boundary normal flux. I guess gradS*n = 0 doesn’t 
>> necessarily mean that the pressure gradient is also zero, but I haven’t 
>> thought about it in detail.
>> 
>>> 
>>>> AFAIK, if you want some mixed Dirichlet/Neumann condition, i.e. fix the 
>>>> pressure + enforce a flux for the second equation, you need to convert the 
>>>> Dirichlet condition to a Robin condition and implement it in “neumann”.
>>> 
>>> Yes, so I was told by DuMux: "Mixed boundary conditions. Use pure boundary 
>>> conditions by converting Dirichlet BCs to Robin BCs". But how exactly do I 
>>> do that? I do not see this kind of thing in the examples (or I do not know 
>>> what to look for). Is it the addRobinFluxDerivatives() call that I should 
>>> implement somehow?
>> 
>> I think you are using cell-centred TPFA. So Dirichlet boundaries are weakly 
>> enforced anyway. You fix (mentally or in the math) the Dirichlet values on 
>> the boundary face and then compute (implement) the normal flux using the 
>> cell-values and the boundary values in the neumann() function of the Problem 
>> class.
>> You get the cell values in the neumann function as “const auto& volVars = 
>> elemVolVars[scvf.insideScvIdx()]” and then e.g. “volVars.pressure(phaseIdx)” 
>> or “volVars.saturation(phaseIdx)”. 
>> I don’t know a 2p example right now, but for example 
>> https://git.iws.uni-stuttgart.de/dumux-repositories/dumux/-/blob/master/test/porousmediumflow/richards/implicit/nonisothermal/evaporation/problem.hh
>>  
>> <https://git.iws.uni-stuttgart.de/dumux-repositories/dumux/-/blob/master/test/porousmediumflow/richards/implicit/nonisothermal/evaporation/problem.hh>
>> implements an energy flux over the boundary like this.
>> 
>> You can implement any boundary flux you want there. But there is 
>> unfortunately no “preimplemented boundary conditions” since the exact 
>> boundary flux highly depends on the considered scenario/application. If you 
>> have suggestion for a “standardized” boundary condition type that is very 
>> useful,
>> feel free to open an issue or suggest an implementation as a merge request: 
>> https://git.iws.uni-stuttgart.de/dumux-repositories/dumux 
>> <https://git.iws.uni-stuttgart.de/dumux-repositories/dumux>. Except of 
>> course Neumann no-flow (which is returning 0.0 in the Problem::neumann 
>> function) or Dirichlet (which is returning the Dirichlet values in the 
>> Problem::dirichlet function).
>> 
>> addRobinFluxDerivatives() is only needed if you want to implement an exact / 
>> analytical Jacobian matrix for your problem. In the default setting DuMux 
>> computes the Jacobian using numerical differentiation and doesn’t need any 
>> user-implemented derivatives.
>> 
>> Hope this helps
>> Timo
>>  
>> 
>>> 
>>> Best regards,
>>> 
>>> Dmitry
>>> 
>>> 
>>> 
>>> On 25.03.2020 22:10, Timo Koch wrote:
>>>> Dear Dmitry,
>>>> 
>>>> depends on what you understand under outflow. The outflow boundary 
>>>> condition is not implemented anymore for the implicit 2p model — exactly 
>>>> for the reason you mention: it’s not really clear what it means.
>>>> I’m wondering why you don’t get an error. I’ll open an issue so that we 
>>>> look at that and that users get a better error message in the future.
>>>> 
>>>> I guess you know better which condition you want in your example, but 
>>>> generally you have two options in Dumux:
>>>> You can set Dirichlet (setAllDirichlet()), depending on your chosen 
>>>> formulation that means setting e.g. pn and Sw., or you can set 
>>>> Neumann/Robin (setAllNeumann).
>>>> For Neumann what we expect from the function “neumann” in the problem 
>>>> class that it returns a boundary flux for both equations (NumEqVector).
>>>> By “flux" I mean the entire term in the divergence. You can look at the 
>>>> documentation for the 2p model: https://dumux.org/doxygen/3.1/a18570.html 
>>>> <https://dumux.org/doxygen/3.1/a18570.html> for the equations implemented.
>>>> 
>>>> AFAIK, if you want some mixed Dirichlet/Neumann condition, i.e. fix the 
>>>> pressure + enforce a flux for the second equation, you need to convert the 
>>>> Dirichlet condition to a Robin condition and implement it in “neumann”.
>>>> Neumann in Dumux also includes “solution-dependent” fluxes (i.e. 
>>>> Robin/Cauchy-type boundary conditions). But depending on your setup, just 
>>>> setting one saturation to 0/1 at the outlet may also give you a good 
>>>> result.
>>>> 
>>>> Best regards,
>>>> Timo
>>>> 
>>>> 
>>>> 
>>>> -- 
>>>> _________________________________________________
>>>> 
>>>> Timo Koch                                      phone: +49 711 685 64676
>>>> IWS, Universität Stuttgart                  fax:   +49 711 685 60430
>>>> Pfaffenwaldring 61         email: [email protected] 
>>>> <mailto:[email protected]>
>>>> D-70569 Stuttgart             url: www.iws.uni-stuttgart.de/en/lh2/ 
>>>> <http://www.iws.uni-stuttgart.de/en/lh2/>
>>>> _________________________________________________
>>>> 
>>>>> On 25. Mar 2020, at 18:39, Dmitry Pavlov <[email protected] 
>>>>> <mailto:[email protected]>> wrote:
>>>>> 
>>>>> Dear Martin,
>>>>> 
>>>>> Thank you! So I was not basing quite on the right model/example for my 
>>>>> task.
>>>>> 
>>>>> I changed the model to dumux/porousmediumflow/2p, and am now having 
>>>>> another problem with the outflow condition: 
>>>>> setOutflow(Indices::saturationIdx). It worked with the sequential model, 
>>>>> however, to be honest, I do not quite understand what it really means. I 
>>>>> would have understood (mathematically/physically) a Neumann condition 
>>>>> that the derivative of the saturation on the outlet boundary must be 
>>>>> zero, but I suspect setOutflow is a different thing, and I never saw a 
>>>>> Neumann condition on saturation in DuMux examples in any form.
>>>>> 
>>>>> Anyway, setOutflow() seemed to work with the sequential model, but with 
>>>>> the implicit model, it immediately results in 
>>>>> Assemble: ... 
>>>>> dumux/dumux/discretization/cellcentered/tpfa/elementvolumevariables.hh:322:
>>>>>  int Dumux::CCTpfaElementVolumeVariables<GVV, false>::getLocalIdx_(int) 
>>>>> const ... Assertion `it != volVarIndices_.end() && "Could not find the 
>>>>> current volume variables for volVarIdx!"' failed.
>>>>> 
>>>>> When I remove setOutflow(), and set e. g. Dirichlet conditions on both 
>>>>> ends, the assertion failure is gone, so I guess that setOutflow() causes 
>>>>> trouble.
>>>>> Could you recommend how do I specify the outflow boundary condition in my 
>>>>> case?
>>>>> 
>>>>> Best regards,
>>>>> 
>>>>> Dmitry
>>>>> 
>>>>> 
>>>>> On 25.03.2020 20:12, Martin Schneider wrote:
>>>>>> Dear Dmitry,
>>>>>> 
>>>>>> but then I don't see the point why you want to use a compositional 
>>>>>> model. 
>>>>>> Couldn't you simply use the dumux/porousmediumflow/2p model?
>>>>>> 
>>>>>> Best regards,
>>>>>> Martin
>>>>>> 
>>>>>> On 25.03.20 18:02, Dmitry Pavlov wrote:
>>>>>>> Dear Martin,
>>>>>>> 
>>>>>>> Thank you for a quick reply. I am currently dealing with a two-phase 
>>>>>>> two-component Darcy flow of water and oil. Both liquids are 
>>>>>>> incompressible and have constant viscosity. Gravity, thermal transfer 
>>>>>>> and diffusion are not considered. Custom functions are provided for 
>>>>>>> capillary pressure, relative permeability for water, and relative 
>>>>>>> permeability for oil. Is more information needed?
>>>>>>> Initially, I wrote my program basing on the example 
>>>>>>> test/porousmediumflow/2p/sequential/test_mpfa2p, and it worked; now I 
>>>>>>> am exploring the implicit method option.
>>>>>>> 
>>>>>>> Best regards,
>>>>>>> Dmitry
>>>>>>> 
>>>>>>> 
>>>>>>> On 25.03.2020 19:50, Martin Schneider wrote:
>>>>>>>> Dear Dmitry,
>>>>>>>> 
>>>>>>>> the compositional models in Dumux take into account diffusion effects 
>>>>>>>> between the different fluid phases,
>>>>>>>> which is why you can't use an immiscible fluidsystem. 
>>>>>>>> 
>>>>>>>> Could you maybe send some more details about your equations you want 
>>>>>>>> to solve
>>>>>>>> and about the phases and components that are involved.
>>>>>>>> 
>>>>>>>> Best regards,
>>>>>>>> Martin
>>>>>>>> 
>>>>>>>> On 25.03.20 17:36, Dmitry Pavlov wrote:
>>>>>>>>> Dear Martin,
>>>>>>>>> 
>>>>>>>>> Following your advice, I decided to try an implicit model of my 
>>>>>>>>> problem, namely a PorousMediumFlowProblem, following the example in 
>>>>>>>>> tests/porousmediumflow/2pnc/implicit/diffusion.
>>>>>>>>> 
>>>>>>>>> I managed to compile it, but am getting a runtime error upon the 
>>>>>>>>> start and, honestly, have no idea what to do about it. I will 
>>>>>>>>> appreciate any pointers.
>>>>>>>>> For convenience, I will describe my task again. I am trying to 
>>>>>>>>> simulate a two-phase two-component flow on a 2D rectangular area with 
>>>>>>>>> one border being the influx and the opposite one being the outflow 
>>>>>>>>> (setOutflow()). Two other two borders are walls, of course, with 
>>>>>>>>> setAllNeumann().
>>>>>>>>> The fluid system is realized as a TwoPImmiscible of two OnePLiquid-s. 
>>>>>>>>> The first of the liquids is the wetting phase, the second is the 
>>>>>>>>> nonwetting phase.
>>>>>>>>> 
>>>>>>>>> Now, the error that I am getting is
>>>>>>>>> 
>>>>>>>>> Dune reported error: Dune::InvalidStateException 
>>>>>>>>> [binaryDiffusionCoefficient:.../dumux/dumux/material/fluidsystems/2pimmiscible.hh:462]:
>>>>>>>>>  Binary diffusion coefficients of components are meaningless if 
>>>>>>>>> immiscibility is assumed
>>>>>>>>> I never any "binary diffusion coefficients" directly in my files, so 
>>>>>>>>> it comes from the method.
>>>>>>>>> 
>>>>>>>>> I have doubts about the following line of my program
>>>>>>>>>     values.setState(Indices::bothPhases);
>>>>>>>>> which is present in dirichletAtPos(), neumannAtPos(), and 
>>>>>>>>> initialAtPos() implementation. I do not really understand what does 
>>>>>>>>> this setState() setting mean, so I kind of guessed.
>>>>>>>>> 
>>>>>>>>> 
>>>>>>>>> Am I on the right path at all here? I can provide more details, of 
>>>>>>>>> course, if that helps to understand what is going on.
>>>>>>>>> 
>>>>>>>>> Best regards,
>>>>>>>>> 
>>>>>>>>> Dmitry
>>>>>>>>> 
>>>>>>>>> 
>>>>>>>>> On 23.03.2020 20:59, Martin Schneider wrote:
>>>>>>>>>> Dear Dmitry,
>>>>>>>>>> 
>>>>>>>>>> most of the Dumux users are using the fully-coupled (implicit) 
>>>>>>>>>> models. 
>>>>>>>>>> The sequential methods have not really been further developed (at 
>>>>>>>>>> least the discretization schemes) 
>>>>>>>>>> over the last few years. 
>>>>>>>>>> 
>>>>>>>>>> However, if you want to use the MPFA-L method you have to use the 
>>>>>>>>>> sequential method.
>>>>>>>>>> The implicit models only support the MPFA-O scheme. 
>>>>>>>>>> 
>>>>>>>>>> Best regards,
>>>>>>>>>> Martin 
>>>>>>>>>> 
>>>>>>>>>> On 23.03.20 18:53, Dmitry Pavlov wrote:
>>>>>>>>>>> Dear Martin,
>>>>>>>>>>> 
>>>>>>>>>>> I guess not, I took the first one that worked, after numerous 
>>>>>>>>>>> failed attempts with other methods taken from random examples, 
>>>>>>>>>>> though the said attempts might have failed because of me being a 
>>>>>>>>>>> non-experienced user, and not because of the nature of the methods. 
>>>>>>>>>>> Would you suggest to use an implicit method instead?
>>>>>>>>>>> 
>>>>>>>>>>> Best regards,
>>>>>>>>>>> 
>>>>>>>>>>> Dmitry
>>>>>>>>>>> On 23.03.2020 20:48, Martin Schneider wrote:
>>>>>>>>>>>> Dear Dmitry,
>>>>>>>>>>>> 
>>>>>>>>>>>> if I remember correctly, the sequential MPFA-L scheme requires 
>>>>>>>>>>>> that each vertex is surrounded by 4 quadrilaterals (in 2D). 
>>>>>>>>>>>> If this is not the case, the construction of interaction volumes 
>>>>>>>>>>>> fails.
>>>>>>>>>>>> 
>>>>>>>>>>>> Is there any reason why you are using the sequential method? 
>>>>>>>>>>>> 
>>>>>>>>>>>> Best regards,
>>>>>>>>>>>> Martin
>>>>>>>>>>>> 
>>>>>>>>>>>> On 23.03.20 18:42, Dmitry Pavlov wrote:
>>>>>>>>>>>>> Dear Kilian,
>>>>>>>>>>>>> 
>>>>>>>>>>>>>> 1.) You should be able to use more recent versions of gmsh if 
>>>>>>>>>>>>>> you go to File -> Export. If you name your grid *.msh
>>>>>>>>>>>>>> 
>>>>>>>>>>>>>> there will be an option to chose the older Version II ASCII 
>>>>>>>>>>>>>> format, readable by Dumux.
>>>>>>>>>>>>> Thank you!
>>>>>>>>>>>>>> 2.) Have you tried an unstructured grid with quadrilaterals? 
>>>>>>>>>>>>>> According to the header where the error comes from,
>>>>>>>>>>>>>> only quadrilaterals are supported (at least, this is what some 
>>>>>>>>>>>>>> comments in the code suggest). You can force gmsh to
>>>>>>>>>>>>>> use only quads (also unstructured ones).
>>>>>>>>>>>>> I just tried an unstructured quad grid and got the same error. 
>>>>>>>>>>>>> However, a gmsh-generated structured quad grid went fine, which 
>>>>>>>>>>>>> means that there is no trouble with gmsh itself.
>>>>>>>>>>>>> 
>>>>>>>>>>>>> Good news: using FVPressureTwoPAdaptive instead of 
>>>>>>>>>>>>> FvMpfaL2dPressureTwoP allowed me to use an gmsh-generated 
>>>>>>>>>>>>> unstructured grid (both triangular and quadrilateral). I am 
>>>>>>>>>>>>> wondering whether MPFA-L in DuMux "officially" requires a 
>>>>>>>>>>>>> structured quad grid. 
>>>>>>>>>>>>> 
>>>>>>>>>>>>> Best regards,
>>>>>>>>>>>>> 
>>>>>>>>>>>>> Dmitry
>>>>>>>>>>>>> 
>>>>>>>>>>>>> 
>>>>>>>>>>>>> 
>>>>>>>>>>>> 
>>>>>>>>>> 
>>>>>>>> 
>>>>>> 
>>>>> _______________________________________________
>>>>> Dumux mailing list
>>>>> [email protected] <mailto:[email protected]>
>>>>> https://listserv.uni-stuttgart.de/mailman/listinfo/dumux 
>>>>> <https://listserv.uni-stuttgart.de/mailman/listinfo/dumux>
>>>> 
>>> _______________________________________________
>>> Dumux mailing list
>>> [email protected] <mailto:[email protected]>
>>> https://listserv.uni-stuttgart.de/mailman/listinfo/dumux 
>>> <https://listserv.uni-stuttgart.de/mailman/listinfo/dumux>
>> 

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

Reply via email to