Hi all,<div><br></div><div>I've created a Matshell multiplication that acts on a vector derived from a DMComposite which contains a DMDA and a DMRedundant of size 1 living on process 0.  It does its job on a single process, but on multiple processes, the VecAXPY (see below) fails with the error:</div>
<div><br></div><div><div><div>[0]PETSC ERROR: --------------------- Error Message ------------------------------------</div><div>[0]PETSC ERROR: Invalid argument!</div><div>[0]PETSC ERROR: Scalar value must be same on all processes, argument # 2!</div>
<div>[0]PETSC ERROR: ------------------------------------------------------------------------</div><div>[0]PETSC ERROR: Petsc Development HG revision: 9cac76b7c349ad267a1605750f02b8ef9d350ff1  HG Date: Tue Aug 21 10:07:43 2012 -0500</div>
<div>[0]PETSC ERROR: See docs/changes/index.html for recent updates.</div><div>[0]PETSC ERROR: See docs/faq.html for hints about trouble shooting.</div><div>[0]PETSC ERROR: See docs/index.html for manual pages.</div><div>
[0]PETSC ERROR: ------------------------------------------------------------------------</div><div>[0]PETSC ERROR: ./ex3 on a arch-linu named lester by kg Tue Aug 21 12:36:10 2012</div><div>[0]PETSC ERROR: Libraries linked from /home/kg/lib/petsc-dev/arch-linux2-c-debug/lib</div>
<div>[0]PETSC ERROR: Configure run at Tue Aug 21 11:54:03 2012</div><div>[0]PETSC ERROR: Configure options PETSC_ARCH=arch-linux2-c-debug</div><div>[0]PETSC ERROR: ------------------------------------------------------------------------</div>
<div>[0]PETSC ERROR: VecAXPY() line 550 in src/vec/vec/interface/rvector.c</div><div>[0]PETSC ERROR: ArcContShellMult() line 1049 in ex3.c</div><div>[0]PETSC ERROR: MatMult() line 2146 in src/mat/interface/matrix.c</div><div>
[0]PETSC ERROR: PCApplyBAorAB() line 621 in src/ksp/pc/interface/precon.c</div><div>[0]PETSC ERROR: KSPGMRESCycle() line 155 in src/ksp/ksp/impls/gmres/gmres.c</div><div>[0]PETSC ERROR: KSPSolve_GMRES() line 231 in src/ksp/ksp/impls/gmres/gmres.c</div>
<div>[0]PETSC ERROR: KSPSolve() line 446 in src/ksp/ksp/interface/itfunc.c</div><div>[0]PETSC ERROR: main() line 425 in ex3.c</div><div>--------------------------------------------------------------------------</div><div>
MPI_ABORT was invoked on rank 0 in communicator MPI_COMM_WORLD</div><div>with errorcode 1.</div></div></div><div><br></div><div>Prodding the Matshell shows that a process that the DMRedundant doesn't "live" on has size 0.  I guess the behaviour I expected is that every process would have space for the redundant entry, but all processes would be synchronized with the process where it lives upon access (or at least on a scatter, broadcast, ... something).  </div>
<div><br></div><div>The examples with DMRedundant are cases where only a single process needs to know the "redundant" part (so why is that called redundant?), and do not demonstrate this behaviour.</div><div><br>
</div><div>Do I have to broadcast the information myself before I can scale a global vector?  Is there something that might work better than DMRedundant?  Any suggestions would be welcome.</div><div><br></div><div>Cheers,</div>
<div>Kevin</div><div><br></div><div>Matshell multiplication routine:</div><div><br></div><div><div>/*                                                                                                                                                                                                                                 </div>
<div>  ArcContShellMult - MatMult overload for linear system of a pseudo-arclength continuation                                                                                                                                         </div>
<div>  [du_post;dlam_post] = [dFdu dFdlam; du0 dlam0]*[du;dlam]                                                                                                                                                                         </div>
<div>                                                                                                                                                                                                                                   </div>
<div>   Input Parameters:                                                                                                                                                                                                               </div>
<div>.  A         - The matshell                                                                                                                                                                                                        </div>
<div>.  x         - the state to operate on                                                                                                                                                                                             </div>
<div>                                                                                                                                                                                                                                   </div>
<div>   Output Parameters:                                                                                                                                                                                                              </div>
<div>.  y         - the result                                                                                                                                                                                                          </div>
<div>*/</div><div>PetscErrorCode ArcContShellMult(Mat A,Vec x,Vec y)</div><div>{</div><div>  PetscInt ierr;</div><div>  void *ctx;</div><div>  ArcContCtx *arc_ctx;</div><div><br></div><div>  Vec delta_U, delta_lam, delta_U_post, delta_lam_post;</div>
<div>  PetscScalar *del_lam_arr,lam_on_entry;</div><div>  PetscFunctionBegin;</div><div>  ierr = MatShellGetContext(A,&ctx); CHKERRQ(ierr);</div><div>  arc_ctx = (ArcContCtx*)(ctx);</div><div><br></div><div>  /*                                                                                                                                                                                                                               </div>
<div>    Get the components of the composite input and output                                                                                                                                                                           </div>
<div>  */</div><div>  ierr = DMCompositeGetAccess(arc_ctx->dmc,x,&delta_U,&delta_lam);CHKERRQ(ierr);</div><div>  ierr = DMCompositeGetAccess(arc_ctx->dmc,y,&delta_U_post,&delta_lam_post);CHKERRQ(ierr);</div>
<div><br></div><div>  /*                                                                                                                                                                                                                               </div>
<div>    Get lambda component as usable array                                                                                                                                                                                           </div>
<div>  */</div><div>  ierr = VecGetArray(delta_lam,&del_lam_arr);CHKERRQ(ierr);</div><div>  lam_on_entry = del_lam_arr[0];</div><div><br></div><div>  /*                                                                                                                                                                                                                               </div>
<div>    Multiply and add each of the components                                                                                                                                                                                        </div>
<div>  */</div><div>  ierr = MatMult(arc_ctx->dFdu,delta_U,delta_U_post);CHKERRQ(ierr);</div><div>  ierr = VecAXPY(delta_U_post,del_lam_arr[0],arc_ctx->dFdlam);CHKERRQ(ierr);</div><div><br></div><div>  ierr = VecDot(arc_ctx->du0,delta_U,&del_lam_arr[0]);CHKERRQ(ierr);</div>
<div>  del_lam_arr[0] += arc_ctx->dlam0*lam_on_entry;</div><div><br></div><div>  /*                                                                                                                                                                                                                               </div>
<div>    Gather components back to the y vector                                                                                                                                                                                         </div>
<div>  */</div><div>  ierr = VecRestoreArray(delta_lam,&del_lam_arr);CHKERRQ(ierr);</div><div>  ierr = VecCopy(delta_lam,delta_lam_post);</div><div>  ierr = DMCompositeRestoreAccess(arc_ctx->dmc,x,&delta_U,&delta_lam);CHKERRQ(ierr);</div>
<div>  ierr = DMCompositeRestoreAccess(arc_ctx->dmc,y,&delta_U_post,&delta_lam_post);CHKERRQ(ierr);</div><div>  PetscFunctionReturn(0);</div><div>}</div><div><br></div></div><div><br></div><div><br></div>