Hi Matt, 

Maybe the nonlinear poro-viscoelastic code Jean-Paul and I contributed to 
the Code Gallery may help you decide whether adding another phase would be 
feasible for your application. Our formulation models a biphasic material 
that consists in a nonlinear (finite-strain) viscoelastic solid and a pore 
fluid. We solve monolithically for solid displacements and fluid pressure 
with a direct solver. The whole code structure is based on that of step-44, 
so it should be relatively easy to navigate if you're already familiar with 
step-44.

Best,
Ester

El dia dimecres, 5 d’octubre de 2022 a les 15:30:49 UTC+2, 
mjri...@gmail.com va escriure:

> This was very helpful! 
>
> I am trying to understand the process so I could potentially extend this 
> formulation for biphasic materials. The way it handles incompressibility is 
> perfect for high water content tissues. My concern revolves around adding 
> another phase and explicitly solving for the fluid displacements in the 
> system. I was going to follow a similar strategy, but I am worried that 
> this approach is very specific. 
>
> IMHO the key feature for this approach is the ability to condense out the 
> pressure and volumetric ratio variables and recover the displacement only 
> formulation but with anti-locking properties baked in as opposed to 
> starting with a pure displacement formulation. I would like to preserve 
> that if possible. I am just worried I will wind up in a hole I cannot climb 
> out of...
>
> These comments were very helpful. One quick question if you just 
> proceed without the condensation would that still work? In one of the 
> tutorials, If I abandoned the notion of using the CG solver and used 
> something like GMRES can I avoid some of the manipulation after I have the 
> linearized system?  
>
> Thanks for the prompt response and thoughtful answers! This stuff is hard 
> and fascinating at the same time. 
>
> Matt 
>
> On Tue, 4 Oct 2022 at 14:03, Jean-Paul Pelteret <jppel...@gmail.com> 
> wrote:
>
>> Hi Matthew,
>>
>> I'm glad that you find step-44 to be a useful tutorial! Let me try to 
>> answer your questions directly.
>>
>> *My first question deals with the statement "The Euler-Lagrange equations 
>> corresponding to the residual"*
>>
>> *Directly above this sentence is the residual, whose derivation I 
>> understand. Where I am lost is that s tact on to equations. I have** 
>> only one residual equation. I cannot bridge that disconnect. *
>>
>> This is just another way of saying that the three equations listed there 
>> (these E-L equations) are the strong form of the governing equations. 
>> Basically, if you take each of these equations (along the way, modifying 
>> the definition of the stress in the equilibrium equation), test them with 
>> the appropriate test function and sum up the three residual contributions 
>> then you recover the (total) residual, or stationary point of the residual, 
>> that is listed above. The point is that it not necessarily so straight 
>> forward to go from the strong form to the weak form for this mixed 
>> formulation, so identifying the conservation equations a-postori is a 
>> helpful sanity check here. They seem to align with what we're trying to do 
>> here.  
>>
>>
>> *I am completely lost here. What is the significance of p and J not 
>> having derivatives on them that makes it "easy" to solve for those terms in 
>> isolation? *
>>
>> Well, that's a valid point. I can't quite recall what exactly we were 
>> trying to identify with this comment. I'll think about it, as that seems to 
>> be a point that we could improve on in the documentation. That there is no 
>> K_{pp} contribution is significant, because it makes condensing out the p 
>> and J fields easier. Maybe we meant to refer to the lack of contribution to 
>> K_{pp} (as there is no second derivative involving a variation and 
>> linearisation of p.
>>
>> *I am also confused how K_pJ, K_Jp and K_JJ  form a block diagonal 
>> matrix. I could get there if I ignore the top row but then the equations 
>> below do not make sense I think. Some detail on this part of the process 
>> would be great. *
>>
>> So this is more easy to explain. We specifically choose discontinuous 
>> shape function to discretise these fields. As there are no interface/flux 
>> contributions, all local element contributions for these terms will remain 
>> local and you therefore end up with an assembled matrix for these 
>> contributions that has a block-like structure. The K_{JJ} matrix is 
>> evidently block-diagonal, as it is a field that couples with itself. As for 
>> the coupling matrices K_{pJ} and K_{Jp}, they are block-diagonal because we 
>> chose *exactly* the same discretisation for both fields (i.e. the shape 
>> functions and polynomial order match).
>>
>> Does that make sense? 
>>
>> Best,
>>
>> Jean-Paul
>> On 2022/09/30 19:16, Matthew Rich wrote:
>>
>> Hi, 
>>
>> Step-44 has a lot going on and really sets you up to tackle a variety of 
>> real world problems. There is a significant jump in complexity between this 
>> tutorial and steps 8,17 & 18. 
>>
>> I have been reading the references and I am about there but need few 
>> clarifications for why things are the way they are. 
>>
>> Any assistance would be much appreciated...
>>
>> My first question deals with the statement "The Euler-Lagrange equations 
>> corresponding to the residual"
>>
>> Directly above this sentence is the residual, whose derivation I 
>> understand. Where I am lost is that s tact on to equations. I have only one 
>> residual equation. I cannot bridge that disconnect. 
>>
>> My other question has to deal with the solving procedure for K (see 
>> snippet below) 
>> [image: Ksolve.png]
>>
>>
>> I am completely lost here. What is the significance of p and J not having 
>> derivatives on them that makes it "easy" to solve for those terms in 
>> isolation? I am also confused how K_pJ, K_Jp and K_JJ  form a block 
>> diagonal matrix. I could get there if I ignore the top row but then the 
>> equations below do not make sense I think. Some detail on this part of the 
>> process would be great. 
>>
>> Thanks in advance, 
>>
>> Matt
>> -- 
>> The deal.II project is located at http://www.dealii.org/
>> For mailing list/forum options, see 
>> https://groups.google.com/d/forum/dealii?hl=en
>> --- 
>> You received this message because you are subscribed to the Google Groups 
>> "deal.II User Group" group.
>> To unsubscribe from this group and stop receiving emails from it, send an 
>> email to dealii+un...@googlegroups.com.
>> To view this discussion on the web visit 
>> https://groups.google.com/d/msgid/dealii/a090a200-fa7a-44d9-ba59-6c5e6064fd5bn%40googlegroups.com
>>  
>> <https://groups.google.com/d/msgid/dealii/a090a200-fa7a-44d9-ba59-6c5e6064fd5bn%40googlegroups.com?utm_medium=email&utm_source=footer>
>> .
>>
>> -- 
>> The deal.II project is located at http://www.dealii.org/
>> For mailing list/forum options, see 
>> https://groups.google.com/d/forum/dealii?hl=en
>> --- 
>>
> You received this message because you are subscribed to a topic in the 
>> Google Groups "deal.II User Group" group.
>> To unsubscribe from this topic, visit 
>> https://groups.google.com/d/topic/dealii/FSaYNTtXNII/unsubscribe.
>> To unsubscribe from this group and all its topics, send an email to 
>> dealii+un...@googlegroups.com.
>> To view this discussion on the web visit 
>> https://groups.google.com/d/msgid/dealii/b5c0a491-7b87-8a56-6c08-36ff5d345377%40gmail.com
>>  
>> <https://groups.google.com/d/msgid/dealii/b5c0a491-7b87-8a56-6c08-36ff5d345377%40gmail.com?utm_medium=email&utm_source=footer>
>> .
>>
>

-- 
The deal.II project is located at http://www.dealii.org/
For mailing list/forum options, see 
https://groups.google.com/d/forum/dealii?hl=en
--- 
You received this message because you are subscribed to the Google Groups 
"deal.II User Group" group.
To unsubscribe from this group and stop receiving emails from it, send an email 
to dealii+unsubscr...@googlegroups.com.
To view this discussion on the web visit 
https://groups.google.com/d/msgid/dealii/c37e6390-9221-4ac0-8459-dc20ccfda317n%40googlegroups.com.

Reply via email to