[petsc-users] Getting a local vector with ghosts knowing ghost indices
Matthew Knepley
knepley at gmail.com
Tue Jan 31 09:26:21 CST 2017
On Sun, Jan 29, 2017 at 3:45 AM, Praveen C <cpraveen at gmail.com> wrote:
> Dear all
>
> In my problem, I know the ghost vertices in each partition since I have
> used metis to partition my unstructured grid. I want to get a vector ul
> with ghosts, from a global vector ug.
>
> I will use ul to assemble my non-linear rhs vector.
>
> Is the following approach optimal for this purpose ?
>
This should work and be scalable
Matt
> PetscInt nvar; // number of variables at each vertex
> PetscInt nvl; // number of local vertices in this partition
> PetscInt nvg; // number of ghost vertices for this partition
> PetscInt *vghost; // global index of ghost vertices for this partition
> PetscReal **u; // size u[nvl+nvg][nvar]
> Vec ul; // vector with ghosts
> Vec ug; // global vector without ghosts
>
> VecCreate(PETSC_COMM_WORLD, &ug);
> VecSetSizes(ug, nvar*nvl, PETSC_DECIDE);
> VecSetFromOptions(ug);
>
> VecCreateGhostBlockWithArray(PETSC_COMM_WORLD, nvar, nvar*nvl,
> PETSC_DECIDE, nvg, vghost,
> &u, ul);
>
> // Following would be inside RHSFunction
> double **array;
> VecGetArray2d(ug, nvl, nvar, 0, 0, &array);
> for(unsigned int i=0; i<nvl; ++i)
> for(unsigned int j=0; j<nvar; ++j)
> u[i][j] = array[i][j];
> VecRestoreArray2d(ug, nvl, nvar, 0, 0, &array);
>
> // Fill ghost values
> VecGhostUpdateBegin(ul, INSERT_VALUES, SCATTER_FORWARD);
> VecGhostUpdateEnd (ul, INSERT_VALUES, SCATTER_FORWARD);
>
> // now use u[][] to compute local part of rhs vector
>
> Thanks
> praveen
>
--
What most experimenters take for granted before they begin their
experiments is infinitely more interesting than any results to which their
experiments lead.
-- Norbert Wiener
-------------- next part --------------
An HTML attachment was scrubbed...
URL: <http://lists.mcs.anl.gov/pipermail/petsc-users/attachments/20170131/c6f27eab/attachment.html>
More information about the petsc-users
mailing list