> On 2. Apr 2020, at 12:40, Dmitry Pavlov <[email protected]> wrote:
> 
> Hello,
> 
> As Timo said,
> 
>> since DuMux 3.0 the “standard" model (what you called implicit method) 
>> actually don’t care if it’s implicit Euler or explicit Euler or sequential.
>> sou can also setup an IMPES scheme in the “new” way of doing things
> 
> So maybe I will just do that. I would like to study how the IMPES scheme will 
> work with my boundary and initial conditions, as compared to the fully 
> implicit one.
> 
> My question is: is there an example of anything resembling a 2p (immiscible, 
> incompressible) problem solved by IMPES method with DuMux 3.0 code?

Hi Dmitry,

I don’t know of an example. For how to do a sequential solver you could look at 
the test/porousmediumflow/tracer/1ptracer examples which first solves a 1p 
problem then a transport problem with the resulting velocity field.
But implementing IMPES will be a bit of a bigger thing. Basic steps would be to 
create two new models, one for the pressure equation, one for the saturation 
equation and solve then after each other in the main file.
Models means implementing volumevariables/localresidual/model.hh just like for 
the 2p model. And of course you need to transfer data between the two models (a 
simple example is 1ptracer).

I can’t promise that there won’t be obstacles on the way or you’ll possibly 
uncover some bugs. (We are highly appreciate bug reports I you found one.)

> 
> 
> I did try some things that I was able to find, namely this line in main.cc:
> 
>     using Assembler = FVAssembler<TypeTag, DiffMethod::numeric>;
> 
> I added "false" to the end to make time discetization explicit
> 
>     using Assembler = FVAssembler<TypeTag, DiffMethod::numeric, false>;
> 
> and got the following translation error:
> 
> .../dumux/dumux/assembly/cclocalassembler.hh:397:58: error: no matching 
> function for call to 
> ‘Dumux::NumericDifferentiation::partialDerivative(Dumux::CCLocalAssembler<TypeTag,
>  Assembler, (Dumux::DiffMethod)0, 
> false>::assembleJacobianAndResidualImpl(Dumux::CCLocalAssembler<TypeTag, 
> Assembler, (Dumux::DiffMethod)0, false>::JacobianMatrix&, 
> Dumux::CCLocalAssembler<TypeTag, Assembler, (Dumux::DiffMethod)0, 
> false>::GridVariables&) ...
> 
> NumericDifferentiation::partialDerivative(evalStorage, elemSol[0][pvIdx], 
> partialDeriv, storageResidual,
> ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~^~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
> eps_(elemSol[0][pvIdx], pvIdx), numDiffMethod);
> ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
> 

You actually ran into a bug here that we recently fixed on master. I just 
backported the fix to the releases/3.1 branch. So if you do a git pull this 
error should be gone.
This is the assembler setting you could use for saturation equation, so yes 
good idea.


> OK, I thought, maybe I will use the DiffMethod::analytic, since I have the 
> MyLocalResidual implementation inherited from 
> TwoPIncompressibleLocalResidual, where addRobinFluxDerivatives is implemented 
> for my neumann()  boundary conditions. The build went fine What I got 
> afterwards was the runtime error
> 
> Solve: M deltax^k = rNewton: Caught exception: "NumericalProblem 
> [solveLinearSystem:/home/dpavlov/DUMUX2/DUMUX/dumux/dumux/nonlinear/newtonsolver.hh:385]:
>  Linear solver did not converge"
> 
> Am I on the right path at all? Should I pursue the analytical derivatives or 
> set up the scheme differently with numeric derivatives?
See above

Best wishes
Timo


-- 
_________________________________________________

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/
_________________________________________________


> 
> 
> Best regards,
> 
> Dmitry
> 
> 
> _______________________________________________
> Dumux mailing list
> [email protected]
> 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