<div dir="ltr"><div class="gmail_extra"><div class="gmail_quote">On Thu, May 28, 2015 at 10:47 AM, Elias Karabelas <span dir="ltr"><<a href="mailto:elias.karabelas@medunigraz.at" target="_blank">elias.karabelas@medunigraz.at</a>></span> wrote:<br><blockquote class="gmail_quote" style="margin:0 0 0 .8ex;border-left:1px #ccc solid;padding-left:1ex">Dear Members,<br>
<br>
I want to solve a Block System arising from the Discretization of a stabilized Finite Element Formulation of the Stokes System.<br>
<br>
I have the following Block Structure<br>
<br>
A -B^T<br>
B C<br>
<br>
The Preconditioner I intend to use is a block preconditioner of the Form<br>
<br>
A -B^T<br>
S<br>
<br>
where S is an approximation of the Schur Complement. For applying the inverse of the schur complement I want to use a Stabilized Least Squares Commutator in the form<br>
<br>
S^-1 = (B diag(Q)^-1 B^T + C_1)^-1 (B diag(Q)^-1 A diag(Q)^-1 B^T + C_2) (B diag(Q)^-1 B^T + C_1)^-1<br>
<br>
where Q is the mass matrix and C_1 and C_2 are some additional stabilization matrices.<br>
<br>
I got from the Manual, that I can use the PCFieldSplit preconditioner for generating the general Block preconditioner as indicated above. And I also found that I can define some arbitrary PC with PCSHELL. My question is, if it is possible to use PCSHELL to define the action of S^-1 as indicated above.<br></blockquote><div><br></div><div>1) Use FieldSplit is the right PC to start with. Make sure you can do something simple like </div><div><br></div><div> A -B^T</div><div> C + B diag(A)^{-1} B^T</div><div><br></div><div>with it before we do the more complicated thing.</div><div><br></div><div>2) You will want to implement a PC for the (1,1) block. You can use a PCSHELL, which is simpler to setup, but</div><div> that means you will have to manually pull out the FieldSplit KSP and set it. If instead you define your own</div><div> PC implementation, its more boilerplate code, but you could specify this PC from the command line without</div><div> any FieldSplit specific code in your application.</div><div><br></div><div>3) Your PC will get two matrices, the MatSchurComplement, and the preconditioning matrix. If you set Q as the</div><div> preconditioning matrix, or really if you set</div><div><br></div><div> A 0</div><div> 0 Q</div><div> </div><div>as the global preconditioning matrix, then the subsolve for (1,1) will get the Schur Complement and Q, and I think</div><div>that is enough to build your Stabilized LSC PC.</div><div><br></div><div>Let me know if this makes sense to you.</div><div><br></div><div> Thanks,</div><div><br></div><div> Matt</div><div><br></div><blockquote class="gmail_quote" style="margin:0 0 0 .8ex;border-left:1px #ccc solid;padding-left:1ex">
Kind Regards<span class="HOEnZb"><font color="#888888"><br>
Elias Karabelas<br>
<br>
-- <br>
Elias Karabelas, Ph.D.<br>
<br>
Medical University of Graz<br>
Institute of Biophysics<br>
Harrachgasse 21/IV<br>
8010 Graz, Austria<br>
<br>
Phone: <a href="tel:%2B43%20316%20380%207759" value="+433163807759" target="_blank">+43 316 380 7759</a><br>
Email: <a href="mailto:elias.karabelas@medunigraz.at" target="_blank">elias.karabelas@medunigraz.at</a><br>
Web : <a href="http://forschung.medunigraz.at/fodok/staff?name=EliasKarabelas" target="_blank">http://forschung.medunigraz.at/fodok/staff?name=EliasKarabelas</a><br>
<br>
</font></span></blockquote></div><br><br clear="all"><div><br></div>-- <br><div class="gmail_signature">What most experimenters take for granted before they begin their experiments is infinitely more interesting than any results to which their experiments lead.<br>-- Norbert Wiener</div>
</div></div>