<html xmlns:o="urn:schemas-microsoft-com:office:office" xmlns:w="urn:schemas-microsoft-com:office:word" xmlns:m="http://schemas.microsoft.com/office/2004/12/omml" xmlns="http://www.w3.org/TR/REC-html40">
<head>
<meta http-equiv="Content-Type" content="text/html; charset=utf-8">
<meta name="Title" content="">
<meta name="Keywords" content="">
<meta name="Generator" content="Microsoft Word 15 (filtered medium)">
<style><!--
/* Font Definitions */
@font-face
        {font-family:"Cambria Math";
        panose-1:2 4 5 3 5 4 6 3 2 4;}
@font-face
        {font-family:Calibri;
        panose-1:2 15 5 2 2 2 4 3 2 4;}
/* Style Definitions */
p.MsoNormal, li.MsoNormal, div.MsoNormal
        {margin:0in;
        margin-bottom:.0001pt;
        font-size:12.0pt;
        font-family:"Times New Roman";}
a:link, span.MsoHyperlink
        {mso-style-priority:99;
        color:blue;
        text-decoration:underline;}
a:visited, span.MsoHyperlinkFollowed
        {mso-style-priority:99;
        color:purple;
        text-decoration:underline;}
span.EmailStyle17
        {mso-style-type:personal-reply;
        font-family:Calibri;
        color:windowtext;}
span.msoIns
        {mso-style-type:export-only;
        mso-style-name:"";
        text-decoration:underline;
        color:teal;}
.MsoChpDefault
        {mso-style-type:export-only;
        font-size:10.0pt;}
@page WordSection1
        {size:8.5in 11.0in;
        margin:1.0in 1.0in 1.0in 1.0in;}
div.WordSection1
        {page:WordSection1;}
--></style>
</head>
<body bgcolor="white" lang="EN-US" link="blue" vlink="purple">
<div class="WordSection1">
<p class="MsoNormal" style="margin-left:.5in">“Steve: are you using the debug version of PETSc that I built?”<o:p></o:p></p>
<p class="MsoNormal"><span style="font-size:11.0pt;font-family:Calibri"><o:p> </o:p></span></p>
<p class="MsoNormal"><span style="font-size:11.0pt;font-family:Calibri">I am.<o:p></o:p></span></p>
<p class="MsoNormal"><span style="font-size:11.0pt;font-family:Calibri"><o:p> </o:p></span></p>
<p class="MsoNormal" style="margin-left:.5in">“Steve: perhaps you could try setting apar = -1.d0, or something and look at your plots to verify that it is not getting updated. Perhaps it is getting updated with 0.0.  And now that I think about it, this assignment
 should put some reasonable value in "apar" so if was not updated it should not be zero.  Is that correct or could a_par be zero coming into this routine?”<o:p></o:p></p>
<p class="MsoNormal"><span style="font-size:11.0pt;font-family:Calibri"><o:p> </o:p></span></p>
<p class="MsoNormal"><span style="font-size:11.0pt;font-family:Calibri">So, Robert and I tested this previously by creating dummy vectors filled with ‘1’, scattering them to petsc, scattering them out of petsc, and then dumping the output. We were still seeing
 just part of the domain, so I don’t think it’s a problem with the array coming in.<o:p></o:p></span></p>
<p class="MsoNormal"><span style="font-size:11.0pt;font-family:Calibri"><o:p> </o:p></span></p>
<p class="MsoNormal"><span style="font-size:11.0pt;font-family:Calibri">I’m reasonably confident that the original vector makes it to petsc just fine when we ‘scatter_from_xgc’, but the whole thing doesn’t make it back when we scatter out of petsc with ‘scatter_from_xgc.’<o:p></o:p></span></p>
<p class="MsoNormal"><span style="font-size:11.0pt;font-family:Calibri"><o:p> </o:p></span></p>
<p class="MsoNormal"><span style="font-size:11.0pt;font-family:Calibri">The troubling bit about this, which caused me to suggest that it was a race condition, was that if I modify the code a little bit I can get the correct result from the scatter. For example,
 we weren’t always checking the returned error flag. If I add in additional error checks, the output shows the entire field. If I do a VecGetArray … VecRestoreArray, I get the whole domain. If I move all the VecScatterEnds and put in an MPI_Barrier before doing
 the final assignment, I get the whole domain.<o:p></o:p></span></p>
<p class="MsoNormal"><span style="font-size:11.0pt;font-family:Calibri"><o:p> </o:p></span></p>
<p class="MsoNormal"><span style="font-size:11.0pt;font-family:Calibri">It is only with the original code that I don’t get the full domain. My hope was that we were making a mistake in syntax and copying the array before it had been flushed from PETSc’s internal
 buffers to the backing array.<o:p></o:p></span></p>
<p class="MsoNormal"><span style="font-size:11.0pt;font-family:Calibri"><o:p> </o:p></span></p>
<p class="MsoNormal"><span style="font-size:11.0pt;font-family:Calibri">Since that’s not the case, it means it’s something more insidious. Which means, as Barry just suggested, I need to fire up valgrind or some other memory checker.<o:p></o:p></span></p>
<p class="MsoNormal"><span style="font-size:11.0pt;font-family:Calibri"><o:p> </o:p></span></p>
<p class="MsoNormal"><span style="font-size:11.0pt;font-family:Calibri">Thanks,<o:p></o:p></span></p>
<p class="MsoNormal"><span style="font-size:11.0pt;font-family:Calibri">--steve<o:p></o:p></span></p>
<p class="MsoNormal"><span style="font-size:11.0pt;font-family:Calibri"><o:p> </o:p></span></p>
<div>
<div>
<p class="MsoNormal">On 8/25/16, 5:18 PM, "Mark Adams" <<a href="mailto:mfadams@lbl.gov">mfadams@lbl.gov</a>> wrote:<o:p></o:p></p>
</div>
</div>
<div>
<p class="MsoNormal"><o:p> </o:p></p>
</div>
<blockquote style="border:none;border-left:solid #B5C4DF 4.5pt;padding:0in 0in 0in 4.0pt;margin-left:3.75pt;margin-right:0in" id="MAC_OUTLOOK_ATTRIBUTION_BLOCKQUOTE">
<div>
<div>
<div>
<p class="MsoNormal">So, this code works on lots of machines and the first scatter does work.
<o:p></o:p></p>
<div>
<p class="MsoNormal"><o:p> </o:p></p>
</div>
<div>
<p class="MsoNormal">And it sounds like there is nothing wrong with this syntax, there is no known race condition here.<o:p></o:p></p>
</div>
<div>
<p class="MsoNormal"><o:p> </o:p></p>
</div>
<div>
<p class="MsoNormal">Steve: are you using the debug version of PETSc that I built?<o:p></o:p></p>
</div>
<div>
<p class="MsoNormal"><o:p> </o:p></p>
</div>
<div>
<p class="MsoNormal">Note, this assignment:<o:p></o:p></p>
</div>
<div>
<p class="MsoNormal"><o:p> </o:p></p>
</div>
<div>
<div>
<p class="MsoNormal"><span style="font-size:9.5pt">  ! scatter solution back - apar<o:p></o:p></span></p>
</div>
<div>
<p class="MsoNormal"><span style="font-size:9.5pt">  apar = a_apar !!!! THIS ONE<o:p></o:p></span></p>
</div>
<div>
<p class="MsoNormal"><span style="font-size:9.5pt">  call VecCreateSeqWithArray(PETSC_COMM_SELF,ione,a_ts%nnode,apar,xxvec(1),ierr);CHKERRQ(ierr)<o:p></o:p></span></p>
</div>
<div>
<p class="MsoNormal"><span style="font-size:9.5pt">  call VecScatterBegin(a_ts%from_petsc(1),a_XX,xxvec(1),INSERT_VALUES,SCATTER_FORWARD,ierr)<o:p></o:p></span></p>
</div>
</div>
<div>
<p class="MsoNormal"><o:p> </o:p></p>
</div>
<div>
<p class="MsoNormal">Is just so that any boundary conditions in nodes that PETSc does not touch, get passed through this routine.<o:p></o:p></p>
</div>
<div>
<p class="MsoNormal"><o:p> </o:p></p>
</div>
<div>
<p class="MsoNormal">Steve: perhaps you could try setting apar = -1.d0, or something and look at your plots to verify that it is not getting updated. Perhaps it is getting updated with 0.0.  And now that I think about it, this assignment should put some reasonable
 value in "apar" so if was not updated it should not be zero.  Is that correct or could a_par be zero coming into this routine?<o:p></o:p></p>
</div>
<div>
<p class="MsoNormal"><o:p> </o:p></p>
</div>
<div>
<p class="MsoNormal">If need be I could add print statements in PETSc to dig into this.<o:p></o:p></p>
</div>
<div>
<p class="MsoNormal"><o:p> </o:p></p>
</div>
</div>
<div>
<p class="MsoNormal"><o:p> </o:p></p>
<div>
<p class="MsoNormal">On Thu, Aug 25, 2016 at 5:04 PM, Mark Adams <<a href="mailto:mfadams@lbl.gov" target="_blank">mfadams@lbl.gov</a>> wrote:<o:p></o:p></p>
<blockquote style="border:none;border-left:solid #CCCCCC 1.0pt;padding:0in 0in 0in 6.0pt;margin-left:4.8pt;margin-right:0in">
<div>
<p class="MsoNormal">Matt is getting confused with Fortran, I just talked with him and he gets it now,<o:p></o:p></p>
</div>
<div>
<div>
<div>
<p class="MsoNormal"><o:p> </o:p></p>
<div>
<p class="MsoNormal">On Thu, Aug 25, 2016 at 4:58 PM, Robert Hager <<a href="mailto:rhager@pppl.gov" target="_blank">rhager@pppl.gov</a>> wrote:<o:p></o:p></p>
<blockquote style="border:none;border-left:solid #CCCCCC 1.0pt;padding:0in 0in 0in 6.0pt;margin-left:4.8pt;margin-right:0in">
<div>
<p class="MsoNormal">Now I am confused.  n1, which is a local variable in this subroutine, is not declared as a pointer or as allocatable here, so there is memory connected to n1. And in my understanding whatever values we get out of XX with the VecScatter
 are written to this memory. After the scatter has finished, we want those values to go into a_n1, which is an in/output of this subroutine. That’s what the assignment a_n1=n1 does.<span style="color:#888888">
<o:p></o:p></span></p>
<div>
<p class="MsoNormal"><span style="color:#888888"><o:p> </o:p></span></p>
</div>
<div>
<p class="MsoNormal"><span style="color:#888888">Robert</span> <o:p></o:p></p>
<div>
<div>
<p class="MsoNormal"><o:p> </o:p></p>
<div>
<p class="MsoNormal"><o:p> </o:p></p>
</div>
<div>
<p class="MsoNormal"><o:p> </o:p></p>
<div>
<blockquote style="margin-top:5.0pt;margin-bottom:5.0pt">
<div>
<p class="MsoNormal">On Aug 25, 2016, at 3:56 PM, Matthew Knepley <<a href="mailto:knepley@gmail.com" target="_blank">knepley@gmail.com</a>> wrote:<o:p></o:p></p>
</div>
<p class="MsoNormal"><o:p> </o:p></p>
<div>
<div>
<div>
<div>
<p class="MsoNormal"><span style="font-size:9.0pt;font-family:Helvetica">On Thu, Aug 25, 2016 at 2:38 PM, Robert Hager <<a href="mailto:rhager@pppl.gov" target="_blank">rhager@pppl.gov</a>> wrote:<o:p></o:p></span></p>
<blockquote style="border:none;border-left:solid #CCCCCC 1.0pt;padding:0in 0in 0in 6.0pt;margin-left:4.8pt;margin-right:0in">
<div>
<p class="MsoNormal"><span style="font-size:9.0pt;font-family:Helvetica">Isn’t xxvec(0) basically a pointer to n1 after the call to VecCreateSeqWithArray? If so, the VecScatter should update the values in n1 and setting a_n1=n1 in the end makes sense.<o:p></o:p></span></p>
</div>
</blockquote>
<div>
<p class="MsoNormal"><span style="font-size:9.0pt;font-family:Helvetica"><o:p> </o:p></span></p>
</div>
<div>
<p class="MsoNormal"><span style="font-size:9.0pt;font-family:Helvetica">xxvec(0) holds a pointer (n1) to the data. That pointer does not change when the scatter writes in to the Vec, so<o:p></o:p></span></p>
</div>
<div>
<p class="MsoNormal"><span style="font-size:9.0pt;font-family:Helvetica">at the end the pointer n1 is still the same, namely a_n1.<o:p></o:p></span></p>
</div>
<div>
<p class="MsoNormal"><span style="font-size:9.0pt;font-family:Helvetica"><o:p> </o:p></span></p>
</div>
<div>
<p class="MsoNormal"><span style="font-size:9.0pt;font-family:Helvetica">   Matt<o:p></o:p></span></p>
</div>
<div>
<p class="MsoNormal"><span style="font-size:9.0pt;font-family:Helvetica"> <o:p></o:p></span></p>
</div>
<blockquote style="border:none;border-left:solid #CCCCCC 1.0pt;padding:0in 0in 0in 6.0pt;margin-left:4.8pt;margin-right:0in">
<div>
<div>
<p class="MsoNormal"><span style="font-size:9.0pt;font-family:Helvetica;color:#888888">Robert</span><span style="font-size:9.0pt;font-family:Helvetica">
<o:p></o:p></span></p>
<div>
<div>
<p class="MsoNormal"><span style="font-size:9.0pt;font-family:Helvetica"><o:p> </o:p></span></p>
<div>
<p class="MsoNormal"><span style="font-size:9.0pt;font-family:Helvetica"><o:p> </o:p></span></p>
</div>
<div>
<p class="MsoNormal"><span style="font-size:9.0pt;font-family:Helvetica"><o:p> </o:p></span></p>
<div>
<blockquote style="margin-top:5.0pt;margin-bottom:5.0pt">
<div>
<p class="MsoNormal"><span style="font-size:9.0pt;font-family:Helvetica">On Aug 25, 2016, at 3:16 PM, Matthew Knepley <<a href="mailto:knepley@gmail.com" target="_blank">knepley@gmail.com</a>> wrote:<o:p></o:p></span></p>
</div>
<p class="MsoNormal"><span style="font-size:9.0pt;font-family:Helvetica"><o:p> </o:p></span></p>
<div>
<div>
<div>
<div>
<p class="MsoNormal"><span style="font-size:9.0pt;font-family:Helvetica">On Thu, Aug 25, 2016 at 1:59 PM, Mark Adams <<a href="mailto:mfadams@lbl.gov" target="_blank">mfadams@lbl.gov</a>> wrote:<o:p></o:p></span></p>
<blockquote style="border:none;border-left:solid #CCCCCC 1.0pt;padding:0in 0in 0in 6.0pt;margin-left:4.8pt;margin-right:0in">
<div>
<div>
<p class="MsoNormal"><span style="font-size:9.0pt;font-family:Helvetica">We have the subroutine below that scatters three vectors.  We have used this code on many machines and it works fine but on one machine data does not get scattered correctly. The first
 scatter looks OK, but it looks like the other two are missing data.<o:p></o:p></span></p>
</div>
<div>
<p class="MsoNormal"><span style="font-size:9.0pt;font-family:Helvetica"><o:p> </o:p></span></p>
</div>
<div>
<p class="MsoNormal"><span style="font-size:9.0pt;font-family:Helvetica">Am I using this correctly?  Should we use VecGetArray in here instead of just using the pointer used for construction? Is there a race condition here that I'm missing?<o:p></o:p></span></p>
</div>
<div>
<p class="MsoNormal"><span style="font-size:9.0pt;font-family:Helvetica"><o:p> </o:p></span></p>
</div>
<div>
<p class="MsoNormal"><span style="font-size:9.0pt;font-family:Helvetica">Thanks,<o:p></o:p></span></p>
</div>
<div>
<p class="MsoNormal"><span style="font-size:9.0pt;font-family:Helvetica">Mark<o:p></o:p></span></p>
</div>
<div>
<p class="MsoNormal"><span style="font-size:9.0pt;font-family:Helvetica"><o:p> </o:p></span></p>
</div>
<div>
<p class="MsoNormal"><span style="font-size:9.0pt;font-family:Helvetica">subroutine scatter_to_xgc(a_ts,a_XX,a_n1,a_apar,a_phi,ierr)<o:p></o:p></span></p>
</div>
<div>
<p class="MsoNormal"><span style="font-size:9.0pt;font-family:Helvetica">  use petscts<o:p></o:p></span></p>
</div>
<div>
<p class="MsoNormal"><span style="font-size:9.0pt;font-family:Helvetica">  use sml_module,only:sml_mype<o:p></o:p></span></p>
</div>
<div>
<p class="MsoNormal"><span style="font-size:9.0pt;font-family:Helvetica">  use xgc_ts_module<o:p></o:p></span></p>
</div>
<div>
<p class="MsoNormal"><span style="font-size:9.0pt;font-family:Helvetica">  implicit none<o:p></o:p></span></p>
</div>
<div>
<p class="MsoNormal"><span style="font-size:9.0pt;font-family:Helvetica">  type(xgc_ts),intent(in)::a_ts<o:p></o:p></span></p>
</div>
<div>
<p class="MsoNormal"><span style="font-size:9.0pt;font-family:Helvetica">  Vec,intent(in)::a_XX<o:p></o:p></span></p>
</div>
<div>
<p class="MsoNormal"><span style="font-size:9.0pt;font-family:Helvetica">  real (kind=8),dimension(a_ts%nnode)::a_n1,a_apar,a_phi<o:p></o:p></span></p>
</div>
<div>
<p class="MsoNormal"><span style="font-size:9.0pt;font-family:Helvetica">  PetscErrorCode,intent(out)::ierr<o:p></o:p></span></p>
</div>
<div>
<p class="MsoNormal"><span style="font-size:9.0pt;font-family:Helvetica">  ! locals<o:p></o:p></span></p>
</div>
<div>
<p class="MsoNormal"><span style="font-size:9.0pt;font-family:Helvetica">  PetscInt,parameter::ione=1<o:p></o:p></span></p>
</div>
<div>
<p class="MsoNormal"><span style="font-size:9.0pt;font-family:Helvetica">  PetscScalar,dimension(a_ts%nnode)::n1,apar,phi<o:p></o:p></span></p>
</div>
<div>
<p class="MsoNormal"><span style="font-size:9.0pt;font-family:Helvetica">  Vec::xxvec(0:2)<o:p></o:p></span></p>
</div>
<div>
<p class="MsoNormal"><span style="font-size:9.0pt;font-family:Helvetica"><o:p> </o:p></span></p>
</div>
<div>
<p class="MsoNormal"><span style="font-size:9.0pt;font-family:Helvetica">  ! scatter solution back - n1<o:p></o:p></span></p>
</div>
<div>
<p class="MsoNormal"><span style="font-size:9.0pt;font-family:Helvetica">  n1 = a_n1<o:p></o:p></span></p>
</div>
</div>
</blockquote>
<div>
<p class="MsoNormal"><span style="font-size:9.0pt;font-family:Helvetica"><o:p> </o:p></span></p>
</div>
<div>
<p class="MsoNormal"><span style="font-size:9.0pt;font-family:Helvetica">Let me be clearer. Here n1 is equal to a_n1, and used to back xxvec(0). It does not change<o:p></o:p></span></p>
</div>
<div>
<p class="MsoNormal"><span style="font-size:9.0pt;font-family:Helvetica">in the course of the function, so why would you have<o:p></o:p></span></p>
</div>
<div>
<p class="MsoNormal"><span style="font-size:9.0pt;font-family:Helvetica"><o:p> </o:p></span></p>
</div>
<div>
<p class="MsoNormal"><span style="font-size:9.0pt;font-family:Helvetica">  a_n1 = n<o:p></o:p></span></p>
</div>
<div>
<p class="MsoNormal"><span style="font-size:9.0pt;font-family:Helvetica"><o:p> </o:p></span></p>
</div>
<div>
<p class="MsoNormal"><span style="font-size:9.0pt;font-family:Helvetica">later on? What am I missing?<o:p></o:p></span></p>
</div>
<div>
<p class="MsoNormal"><span style="font-size:9.0pt;font-family:Helvetica"><o:p> </o:p></span></p>
</div>
<div>
<p class="MsoNormal"><span style="font-size:9.0pt;font-family:Helvetica">   Matt<o:p></o:p></span></p>
</div>
<div>
<p class="MsoNormal"><span style="font-size:9.0pt;font-family:Helvetica"> <o:p></o:p></span></p>
</div>
<blockquote style="border:none;border-left:solid #CCCCCC 1.0pt;padding:0in 0in 0in 6.0pt;margin-left:4.8pt;margin-right:0in">
<div>
<div>
<p class="MsoNormal"><span style="font-size:9.0pt;font-family:Helvetica">  call VecCreateSeqWithArray(PETSC_COMM_SELF,ione,a_ts%nnode,n1,xxvec(0),ierr);CHKERRQ(ierr)<o:p></o:p></span></p>
</div>
<div>
<p class="MsoNormal"><span style="font-size:9.0pt;font-family:Helvetica">  call VecScatterBegin(a_ts%from_petsc(0),a_XX,xxvec(0),INSERT_VALUES,SCATTER_FORWARD,ierr)<o:p></o:p></span></p>
</div>
<div>
<p class="MsoNormal"><span style="font-size:9.0pt;font-family:Helvetica"><o:p> </o:p></span></p>
</div>
<div>
<p class="MsoNormal"><span style="font-size:9.0pt;font-family:Helvetica">  ! scatter solution back - apar<o:p></o:p></span></p>
</div>
<div>
<p class="MsoNormal"><span style="font-size:9.0pt;font-family:Helvetica">  apar = a_apar<o:p></o:p></span></p>
</div>
<div>
<p class="MsoNormal"><span style="font-size:9.0pt;font-family:Helvetica">  call VecCreateSeqWithArray(PETSC_COMM_SELF,ione,a_ts%nnode,apar,xxvec(1),ierr);CHKERRQ(ierr)<o:p></o:p></span></p>
</div>
<div>
<p class="MsoNormal"><span style="font-size:9.0pt;font-family:Helvetica">  call VecScatterBegin(a_ts%from_petsc(1),a_XX,xxvec(1),INSERT_VALUES,SCATTER_FORWARD,ierr)<o:p></o:p></span></p>
</div>
<div>
<p class="MsoNormal"><span style="font-size:9.0pt;font-family:Helvetica"><o:p> </o:p></span></p>
</div>
<div>
<p class="MsoNormal"><span style="font-size:9.0pt;font-family:Helvetica">  ! scatter solution back - phi<o:p></o:p></span></p>
</div>
<div>
<p class="MsoNormal"><span style="font-size:9.0pt;font-family:Helvetica">  phi = a_phi<o:p></o:p></span></p>
</div>
<div>
<p class="MsoNormal"><span style="font-size:9.0pt;font-family:Helvetica">  call VecCreateSeqWithArray(PETSC_COMM_SELF,ione,a_ts%nnode,phi,xxvec(2),ierr);CHKERRQ(ierr)<o:p></o:p></span></p>
</div>
<div>
<p class="MsoNormal"><span style="font-size:9.0pt;font-family:Helvetica">  call VecScatterBegin(a_ts%from_petsc(2),a_XX,xxvec(2),INSERT_VALUES,SCATTER_FORWARD,ierr)<o:p></o:p></span></p>
</div>
<div>
<p class="MsoNormal"><span style="font-size:9.0pt;font-family:Helvetica"><o:p> </o:p></span></p>
</div>
<div>
<p class="MsoNormal"><span style="font-size:9.0pt;font-family:Helvetica">  ! end<o:p></o:p></span></p>
</div>
<div>
<p class="MsoNormal"><span style="font-size:9.0pt;font-family:Helvetica">  call VecScatterEnd(  a_ts%from_petsc(0),a_XX,xxvec(0),INSERT_VALUES,SCATTER_FORWARD,ierr)<o:p></o:p></span></p>
</div>
<div>
<p class="MsoNormal"><span style="font-size:9.0pt;font-family:Helvetica">  a_n1 = n1<o:p></o:p></span></p>
</div>
<div>
<p class="MsoNormal"><span style="font-size:9.0pt;font-family:Helvetica">  call VecDestroy(xxvec(0),ierr)<o:p></o:p></span></p>
</div>
<div>
<p class="MsoNormal"><span style="font-size:9.0pt;font-family:Helvetica"><o:p> </o:p></span></p>
</div>
<div>
<p class="MsoNormal"><span style="font-size:9.0pt;font-family:Helvetica">  call VecScatterEnd(  a_ts%from_petsc(1),a_XX,xxvec(1),INSERT_VALUES,SCATTER_FORWARD,ierr)<o:p></o:p></span></p>
</div>
<div>
<p class="MsoNormal"><span style="font-size:9.0pt;font-family:Helvetica">  a_apar = apar<o:p></o:p></span></p>
</div>
<div>
<p class="MsoNormal"><span style="font-size:9.0pt;font-family:Helvetica">  call VecDestroy(xxvec(1),ierr)<o:p></o:p></span></p>
</div>
<div>
<p class="MsoNormal"><span style="font-size:9.0pt;font-family:Helvetica"><o:p> </o:p></span></p>
</div>
<div>
<p class="MsoNormal"><span style="font-size:9.0pt;font-family:Helvetica">  call VecScatterEnd(  a_ts%from_petsc(2),a_XX,xxvec(2),INSERT_VALUES,SCATTER_FORWARD,ierr)<o:p></o:p></span></p>
</div>
<div>
<p class="MsoNormal"><span style="font-size:9.0pt;font-family:Helvetica">  a_phi = phi<o:p></o:p></span></p>
</div>
<div>
<p class="MsoNormal"><span style="font-size:9.0pt;font-family:Helvetica">  call VecDestroy(xxvec(2),ierr)<o:p></o:p></span></p>
</div>
<div>
<p class="MsoNormal"><span style="font-size:9.0pt;font-family:Helvetica"><o:p> </o:p></span></p>
</div>
<div>
<p class="MsoNormal"><span style="font-size:9.0pt;font-family:Helvetica">  return<o:p></o:p></span></p>
</div>
<div>
<p class="MsoNormal"><span style="font-size:9.0pt;font-family:Helvetica">end subroutine scatter_to_xgc<o:p></o:p></span></p>
</div>
</div>
</blockquote>
</div>
<p class="MsoNormal"><span style="font-size:9.0pt;font-family:Helvetica"><br>
<br clear="all">
<o:p></o:p></span></p>
<div>
<p class="MsoNormal"><span style="font-size:9.0pt;font-family:Helvetica"><o:p> </o:p></span></p>
</div>
<p class="MsoNormal"><span style="font-size:9.0pt;font-family:Helvetica">-- <o:p></o:p></span></p>
<div>
<p class="MsoNormal"><span style="font-size:9.0pt;font-family:Helvetica">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<o:p></o:p></span></p>
</div>
</div>
</div>
</div>
</blockquote>
</div>
<p class="MsoNormal"><span style="font-size:9.0pt;font-family:Helvetica"><o:p> </o:p></span></p>
</div>
</div>
</div>
</div>
</div>
</blockquote>
</div>
<p class="MsoNormal"><span style="font-size:9.0pt;font-family:Helvetica"><br>
<br clear="all">
<o:p></o:p></span></p>
<div>
<p class="MsoNormal"><span style="font-size:9.0pt;font-family:Helvetica"><o:p> </o:p></span></p>
</div>
<p class="MsoNormal"><span style="font-size:9.0pt;font-family:Helvetica">-- <o:p></o:p></span></p>
<div>
<p class="MsoNormal"><span style="font-size:9.0pt;font-family:Helvetica">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<o:p></o:p></span></p>
</div>
</div>
</div>
</div>
</blockquote>
</div>
<p class="MsoNormal"><o:p> </o:p></p>
</div>
</div>
</div>
</div>
</div>
</blockquote>
</div>
<p class="MsoNormal"><o:p> </o:p></p>
</div>
</div>
</div>
</blockquote>
</div>
<p class="MsoNormal"><o:p> </o:p></p>
</div>
</div>
</div>
</blockquote>
</div>
</body>
</html>