Matt's method seems to work well, though instead of editing the actual function 
I put the relevant parts directly into my code. I made the small example 
attached here.


I might look into Star Forests at some point, though it's not really touched 
upon in the manual (I will probably take a look at your paper, 
https://arxiv.org/abs/2102.13018).


Med venlig hilsen / Best regards

Peder

________________________________
Fra: Junchao Zhang <[email protected]>
Sendt: 1. juli 2021 16:38:29
Til: Jed Brown
Cc: Peder Jørgensgaard Olesen; [email protected]
Emne: Re: Sv: [petsc-users] Scatter parallel Vec to sequential Vec on 
non-zeroth process

Peder,
  PETSCSF_PATTERN_ALLTOALL only supports MPI_Alltoall (not Alltoallv), and is 
only used by petsc internally at few places.
  I suggest you can go with Matt's approach. After it solves your problem, you 
can distill an example to demo the communication pattern. Then we can see how 
to efficiently support that in petsc.

  Thanks.
--Junchao Zhang


On Thu, Jul 1, 2021 at 7:42 AM Jed Brown 
<[email protected]<mailto:[email protected]>> wrote:
Peder Jørgensgaard Olesen <[email protected]<mailto:[email protected]>> writes:

> Each process is assigned an indexed subset of the tasks (the tasks are of 
> constant size), and, for each task index, the relevant data is scattered as a 
> SEQVEC to the process (this is done for all processes in each step, using an 
> adaption of the code in Matt's link). This way each process only receives 
> just the data it needs to complete the task. While I'm currently working with 
> very moderate size data sets I'll eventually need to handle something rather 
> more massive, so I want to economize memory where possible and give each 
> process only the data it needs.

>From the sounds of it, this pattern ultimately boils down to MPI_Gather being 
>called P times where P is the size of the communicator. This will work okay 
>when P is small, but it's much less efficient than calling MPI_Alltoall (or 
>MPI_Alltoallv), which you can do by creating one PetscSF that ships the needed 
>data to each task and PETSCSF_PATTERN_ALLTOALL. You can see an example.

https://gitlab.com/petsc/petsc/-/blob/main/src/vec/is/sf/tests/ex3.c#L93-151
#include <stdio.h>
#include <petscvec.h>

int main(int argc, char **argv){
  PetscErrorCode ierr;
  PetscInt n = 5, loc_tgt_size;
  PetscMPIInt rank, size;
  PetscReal vnorm;
  Vec src_vec, tgt_vec;
  PetscRandom rctx;
  IS is;
  VecScatter sctx;

  ierr = PetscInitialize(&argc, &argv, NULL, NULL); CHKERRQ(ierr);

  MPI_Comm_size(PETSC_COMM_WORLD, &size);
  MPI_Comm_rank(PETSC_COMM_WORLD, &rank);

  // Create parallel source vector with random entries
  ierr = VecCreate(PETSC_COMM_WORLD, &src_vec); CHKERRQ(ierr);
  ierr = VecSetType(src_vec, VECSTANDARD); CHKERRQ(ierr);
  ierr = VecSetSizes(src_vec, PETSC_DECIDE, n); CHKERRQ(ierr);
  ierr = PetscRandomCreate(PETSC_COMM_WORLD, &rctx); CHKERRQ(ierr);
  ierr = VecSetRandom(src_vec, rctx); CHKERRQ(ierr);
  ierr = PetscRandomDestroy(&rctx); CHKERRQ(ierr);

  // Create sequential target vector w. zero length on all but one process, zero entries.
  loc_tgt_size = 0;
  if (rank == (size-1)){
    loc_tgt_size = n;
  }
  ierr = VecCreateSeq(PETSC_COMM_SELF, loc_tgt_size, &tgt_vec); CHKERRQ(ierr);
  ierr = VecZeroEntries(tgt_vec); CHKERRQ(ierr);

  // Scatter source vector to target vector on one process
  ierr = ISCreateStride(PETSC_COMM_SELF, loc_tgt_size, 0, 1, &is); CHKERRQ(ierr);
  ierr = VecScatterCreate(src_vec, is, tgt_vec, is, &sctx); CHKERRQ(ierr);
  ierr = VecScatterBegin(sctx, src_vec, tgt_vec, INSERT_VALUES, SCATTER_FORWARD); CHKERRQ(ierr);
  ierr = VecScatterEnd(sctx, src_vec, tgt_vec, INSERT_VALUES, SCATTER_FORWARD); CHKERRQ(ierr);

  // Get and print norm of target vector on each process
  ierr = VecNorm(tgt_vec, NORM_2, &vnorm); CHKERRQ(ierr);
  ierr = PetscPrintf(PETSC_COMM_SELF, "[%i] |v| = %e\n", rank, vnorm); CHKERRQ(ierr);

  ierr = VecDestroy(&src_vec); CHKERRQ(ierr);
  ierr = VecDestroy(&tgt_vec); CHKERRQ(ierr);
  ierr = ISDestroy(&is); CHKERRQ(ierr);
  ierr = VecScatterDestroy(&sctx); CHKERRQ(ierr);

  ierr = PetscFinalize(); CHKERRQ(ierr);

  return(ierr);
}

Reply via email to