<html>
<head>
<meta http-equiv="Content-Type" content="text/html; charset=UTF-8">
</head>
<body text="#000000" bgcolor="#FFFFFF">
<div class="moz-cite-prefix">On 2/19/19 7:40 AM, Matthew Knepley
wrote:<br>
</div>
<blockquote type="cite"
cite="mid:CAMYG4GmbhWj_Q9fFLf6ADY6ZSWvQJLyXKs5qTfB7h43K2T-X2A@mail.gmail.com">
<meta http-equiv="content-type" content="text/html; charset=UTF-8">
<div dir="ltr">
<div dir="ltr">On Tue, Feb 19, 2019 at 1:13 AM Jordan Wagner via
petsc-users <<a href="mailto:petsc-users@mcs.anl.gov"
moz-do-not-send="true">petsc-users@mcs.anl.gov</a>>
wrote:<br>
</div>
<div class="gmail_quote">
<blockquote class="gmail_quote" style="margin:0px 0px 0px
0.8ex;border-left:1px solid
rgb(204,204,204);padding-left:1ex">
<div bgcolor="#FFFFFF">
<p>Hi, <br>
</p>
<p>Over the past few months, I have implemented
dmplex/section in my preexisting finite element code. We
already had our own implementations for things like the
finite element space, boundary conditions, projections,
etc, so I am mostly using dmplex/section at its lowest
level; no DT/DS stuff for the moment. Everything is
working fine at this level, and I am overall very happy
with how it ties into the code that we already have. <br>
</p>
<p>I am currently doing my own preallocation for the
system Jacobian but am wanting to use DMCreateMatrix()
instead, as it seems pointless to keep using my own
crappy function when this one is just sitting there
waiting. However, I dislike how I am currently
implementing this and was hoping to get some suggestions
or some advice. My problem is handling nicely the
Dirichlet dofs. My code takes the approach of first
assembling the entire Jacobian matrix and load vector as
if no constraints are imposed. We then have a function
that applies the essential conditions after assembly and
extracts a submatrix and subvector similar to the way <a
href="https://www.mcs.anl.gov/petsc/petsc-current/src/ksp/ksp/examples/tutorials/ex49.c.html"
target="_blank" moz-do-not-send="true">ex49</a> is
doing. <br>
</p>
<p>Since DMCreateMatrix() automatically partitions out the
constraints specified in the section, I find myself
having to create two nearly identical sections: one that
has constraints and one that does not, and then I clone
the DM and set the default section of each respectively.
I can then use DMCreateMatrix() on the one dm/section
with no constraints to preallocate the larger matrix
that I am assembling. Then I use the dm/section with
constraints to do all of my boundary condition and
projection operations on the previously created matrix
and vector. Essentially, I am using an unconstrained
section only for allocation and a constrained section
for most everything else.<br>
</p>
<p>This process actually works fine for me, but its pretty
ugly, and I am sure there is a better way. I am
wondering if I am missing something that could keep me
from having to go through this process of creating two
sections that differ only in the constraints. It seems
to me if there were an option for
DMCreateMatrix()/DMCreateVector() to either include or
not include the dofs, that would completely solve my
problem. That way, I can use the same section for both
creation and allocation.</p>
<p>So, the main question I have is just whether there is a
function or option that I am missing like DMCreateXXXX()
but with the ability to retain the constrained rows and
columns rather than compressing them out.</p>
</div>
</blockquote>
<div>I do not completely understand the above, but here is how
Section works. A _local_ Section contains all dofs and is
intended for assembly using local Vecs. To communicate with
the solver, you make a _global_ Section using
PetscSectionCreateGlobalSection(). This has an option to
include or exclude constraints. For example, I make a global
Section with constraints when I output for visualization.
Maybe you can use that to do what you want.</div>
<div><br>
</div>
</div>
</div>
</blockquote>
<p>Thanks for the prompt reply and sorry for the unclear initial
question. I think this helps in my understanding of how this tool
was intended to be used. <br>
</p>
<blockquote type="cite"
cite="mid:CAMYG4GmbhWj_Q9fFLf6ADY6ZSWvQJLyXKs5qTfB7h43K2T-X2A@mail.gmail.com">
<div dir="ltr">
<div class="gmail_quote">
<div>What I do not understand is why you do not just use a
global Section without constraints to assemble your matrix,
rather than eliminating the rows/columns later. The global
Section gives back negative indices for unknown that are
constrained or not owned. These are ignored by
MatSetValues() and VecSetValues(). This is how the default
Plex assembly makes operators without constrained variables.</div>
<div><br>
</div>
</div>
</div>
</blockquote>
<p>I have rarely used a global section directly because I have never
been able to make sense of exactly what it was giving me. For
instance as you mentioned (and I have seen written elsewhere), the
global section should give negative indices for unowned or
constrained dofs. This seems to work fine for me with the unowned
points---I clearly get a negative size and offset from the section
on these points. But my problem is for points that have
constraints. I don't see where I would get negative indices from.
The sizes and offsets for these points are not negative. Of
course, from the constraint size and constraint indices on these
points, I can manually find out which should be negative and go
from there. However, that leads me to the next problem I am having
with global section...</p>
<p>For some reason, when I am creating the global section from the
local one, my constrained component indices are getting messed up.
The sizes are right but the order seems to be random and/or
nonsensical. For instance, here is a snippet from my code being
run in serial<br>
</p>
<p>Local section:</p>
<p> (17615) dim 2 offset 794 constrained 0<br>
(17616) dim 2 offset 796 constrained 0 1<br>
(17617) dim 2 offset 798 constrained 0 1<br>
(17618) dim 2 offset 800 constrained 0 1<br>
(17619) dim 2 offset 802 constrained 0 1<br>
(17620) dim 2 offset 804 constrained 0 1<br>
(17621) dim 2 offset 806 constrained 0 1<br>
<br>
</p>
<p>Global section:</p>
<p> (17615) dim 2 offset 80 constrained 1<br>
(17616) dim 2 offset 81 constrained 0 0<br>
(17617) dim 2 offset 81 constrained 0 0<br>
(17618) dim 2 offset 81 constrained 0 0<br>
(17619) dim 2 offset 81 constrained 0 0<br>
(17620) dim 2 offset 81 constrained 0 0<br>
(17621) dim 2 offset 81 constrained 0 1<br>
<br>
</p>
<p>So using the global section, I am unable to even manually tell
which components are constrained. Is this normal for the global
section or is it something going wrong in its creation? or am I
completely misunderstanding something? <br>
</p>
<p>Thanks very much.<br>
</p>
<p><br>
</p>
<blockquote type="cite"
cite="mid:CAMYG4GmbhWj_Q9fFLf6ADY6ZSWvQJLyXKs5qTfB7h43K2T-X2A@mail.gmail.com">
<div dir="ltr">
<div class="gmail_quote">
<div> Thanks,</div>
<div><br>
</div>
<div> Matt</div>
<blockquote class="gmail_quote" style="margin:0px 0px 0px
0.8ex;border-left:1px solid
rgb(204,204,204);padding-left:1ex">
<div bgcolor="#FFFFFF">
<p>Thanks a lot, as always.<br>
</p>
</div>
</blockquote>
</div>
<br clear="all">
<div><br>
</div>
-- <br>
<div dir="ltr" class="gmail_signature">
<div dir="ltr">
<div>
<div dir="ltr">
<div>
<div dir="ltr">
<div>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><br>
</div>
<div><a href="http://www.cse.buffalo.edu/~knepley/"
target="_blank" moz-do-not-send="true">https://www.cse.buffalo.edu/~knepley/</a><br>
</div>
</div>
</div>
</div>
</div>
</div>
</div>
</div>
</blockquote>
</body>
</html>