<div dir="ltr"><div class="gmail_extra"><div class="gmail_quote">On Wed, Aug 6, 2014 at 10:47 PM, Justin Chang <span dir="ltr"><<a href="mailto:jychang48@gmail.com" target="_blank">jychang48@gmail.com</a>></span> wrote:<br>
<blockquote class="gmail_quote" style="margin:0px 0px 0px 0.8ex;border-left-width:1px;border-left-color:rgb(204,204,204);border-left-style:solid;padding-left:1ex"><div dir="ltr"><div><div><div>Thank you for the answer, i have a few more questions:<br>
<br></div>1) In ex62 regarding the functions f0/f1, what's the difference between them? Do f0 correspond with test functions and f1 with the derivatives of the test functions?<br></div></div></div></blockquote><div><br>
</div><div>The residual F is calculated using</div><div><br></div><div>    F(u) = \int_\Omega phi f_0(u, grad u) + grad phi . f_1(u, grad u)</div><div><br></div><div>This formulation is in my slides and also this paper <a href="http://arxiv.org/abs/1309.1204">http://arxiv.org/abs/1309.1204</a>.</div>
<div> </div><blockquote class="gmail_quote" style="margin:0px 0px 0px 0.8ex;border-left-width:1px;border-left-color:rgb(204,204,204);border-left-style:solid;padding-left:1ex"><div dir="ltr"><div><div></div>2) For darcy's law, the velocity boundary conditions for each surface is simply the normal. i.e. u dot n over Gamma_u = <some value>. In DMPlexAddBoundary it seems to me that it constraints every dof corresponding to the mesh point (and field), but is there a way where given a field (say in this case velocity) I only want to prescribe the x velocity on the left and right faces and y velocity on the top and bottom?<br>
</div></div></blockquote><div><br></div><div>You are describing a Neumann condition. If you look at ex12, I show both Dirichlet and Neumann conditions. The Neumann conditions end up</div><div>integrating your <some value> function over the boundary (that is what comes out of integration by parts). So you do that AddBoundary() setting</div>
<div>isEssential to PETSC_FALSE.</div><div> </div><blockquote class="gmail_quote" style="margin:0px 0px 0px 0.8ex;border-left-width:1px;border-left-color:rgb(204,204,204);border-left-style:solid;padding-left:1ex"><div dir="ltr">
<div></div>3) Is there a way/function for DMPlex to output multiple fields to VTK files? Or would the field solutions have to be manually separated out from the original global solution?</div></blockquote><div><br></div><div>
VTK is very limited. I would switch to using HDF5. The fields are automatically split up and there is a Python</div><div>script for generating Xdmf files which make the HDF5 usable from Paraview just like VTK.</div><div><br>
</div><div>If you do not like that, you can try Jed's VTU but I have not tried it yet.</div><div><br></div><div>  Thanks,</div><div><br></div><div>     Matt</div><div> <br></div><blockquote class="gmail_quote" style="margin:0px 0px 0px 0.8ex;border-left-width:1px;border-left-color:rgb(204,204,204);border-left-style:solid;padding-left:1ex">
<div dir="ltr"><div><div><div><div><div><div>Thanks,<br>Justin<br></div></div></div></div></div></div></div><div class="gmail_extra"><br><br><div class="gmail_quote">On Wed, Aug 6, 2014 at 6:27 PM, Matthew Knepley <span dir="ltr"><<a href="mailto:knepley@gmail.com" target="_blank">knepley@gmail.com</a>></span> wrote:<br>

<blockquote class="gmail_quote" style="margin:0px 0px 0px 0.8ex;border-left-width:1px;border-left-color:rgb(204,204,204);border-left-style:solid;padding-left:1ex"><div dir="ltr"><div class="gmail_extra"><div class="gmail_quote">
<div><div>On Wed, Aug 6, 2014 at 5:55 PM, Justin Chang <span dir="ltr"><<a href="mailto:jychang48@gmail.com" target="_blank">jychang48@gmail.com</a>></span> wrote:<br>

<blockquote class="gmail_quote" style="margin:0px 0px 0px 0.8ex;border-left-width:1px;border-left-color:rgb(204,204,204);border-left-style:solid;padding-left:1ex"><div dir="ltr"><div><div><div>Hi,<br><br>So in ex62 of SNES it solves the stokes equation with P2-P1 elements. I want to solve the darcy equation with P2-P1 elements in a similar fashion.<br>


<br></div>The difference between stokes and darcy is that in the weak form of darcy, the velocity-velocity term is v*u whereas in stokes it's grad[v]*grad[u]<br>
<br></div>In ex62, the function g3_uu (lines 172-183) formulates the grad[v]*grad[u] term but how would I do it for Darcy's equation where its simply v*u (assuming drag coefficient is scalar quantity of 1.0)? <br></div>


</div></blockquote><div><br></div></div></div><div>So you would want to define g0 instead of g3, since its the u-u term:</div><div><br></div><div><div>void g0_uu(const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], const PetscReal x[], PetscScalar g3[])</div>


<div>{</div><div>  const PetscInt Ncomp = spatialDim;<br></div><div>  PetscInt       compI;</div><div><br></div><div>  for (compI = 0; compI < Ncomp; ++compI) g0[compI] = u[compI];</div><div>}</div></div><div><br></div>


<div>and then register it in SetupProblem:</div><div><br></div><div><div>  ierr = PetscDSSetJacobian(prob, 0, 0, g0_uu, NULL,  NULL,  NULL);CHKERRQ(ierr);</div></div><div><br></div><div>and do not forget to change the residual definition in f0/f1.</div>


<div><br></div><div>  Thanks,</div><div><br></div><div>     Matt</div><div> </div><blockquote class="gmail_quote" style="margin:0px 0px 0px 0.8ex;border-left-width:1px;border-left-color:rgb(204,204,204);border-left-style:solid;padding-left:1ex">


<div dir="ltr"><div>
</div>Thanks,<br>Justin<span><font color="#888888"><br></font></span></div><span><font color="#888888">
</font></span></blockquote></div><span><font color="#888888"><br><br clear="all"><span class=""><font color="#888888"><div><br></div>-- <br>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
</font></span></font></span></div></div>
</blockquote></div><br></div>
</blockquote></div><br><br clear="all"><div><br></div>-- <br>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>