xref: /petsc/src/vec/is/sf/tests/ex21.c (revision 54a2460711452ea1d7d048b8e4cf631dfb94366e)
1*54a24607SJunchao Zhang static const char help[] = "Test VecScatterCopy() on an SF with duplicated leaves \n\n";
2*54a24607SJunchao Zhang 
3*54a24607SJunchao Zhang #include <petscvec.h>
4*54a24607SJunchao Zhang #include <petscsf.h>
5*54a24607SJunchao Zhang 
6*54a24607SJunchao Zhang /*
7*54a24607SJunchao Zhang   Contributed-by: "Hammond, Glenn E" <glenn.hammond@pnnl.gov>
8*54a24607SJunchao Zhang */
9*54a24607SJunchao Zhang int main(int argc,char* argv[])
10*54a24607SJunchao Zhang {
11*54a24607SJunchao Zhang   PetscMPIInt size;
12*54a24607SJunchao Zhang   PetscInt    n;
13*54a24607SJunchao Zhang   PetscInt    *indices;
14*54a24607SJunchao Zhang   Vec         vec;
15*54a24607SJunchao Zhang   Vec         vec2;
16*54a24607SJunchao Zhang   IS          is;
17*54a24607SJunchao Zhang   IS          is2;
18*54a24607SJunchao Zhang   VecScatter  scatter;
19*54a24607SJunchao Zhang   VecScatter  scatter2;
20*54a24607SJunchao Zhang 
21*54a24607SJunchao Zhang   PetscCall(PetscInitialize(&argc,&argv,NULL,help));
22*54a24607SJunchao Zhang   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
23*54a24607SJunchao Zhang   PetscCheck(size == 1,PETSC_COMM_WORLD,PETSC_ERR_WRONG_MPI_SIZE,"This example requires 1 process");
24*54a24607SJunchao Zhang 
25*54a24607SJunchao Zhang   n = 4;
26*54a24607SJunchao Zhang   PetscCall(PetscMalloc1(n,&indices));
27*54a24607SJunchao Zhang   indices[0] = 0;
28*54a24607SJunchao Zhang   indices[1] = 1;
29*54a24607SJunchao Zhang   indices[2] = 2;
30*54a24607SJunchao Zhang   indices[3] = 3;
31*54a24607SJunchao Zhang   PetscCall(ISCreateGeneral(PETSC_COMM_WORLD,n,indices,PETSC_COPY_VALUES,&is));
32*54a24607SJunchao Zhang   PetscCall(PetscFree(indices));
33*54a24607SJunchao Zhang   PetscCall(VecCreateMPI(PETSC_COMM_WORLD,n,n,&vec));
34*54a24607SJunchao Zhang 
35*54a24607SJunchao Zhang   n = 4;
36*54a24607SJunchao Zhang   PetscCall(PetscMalloc1(n,&indices));
37*54a24607SJunchao Zhang   indices[0] = 0;
38*54a24607SJunchao Zhang   indices[1] = 0;
39*54a24607SJunchao Zhang   indices[2] = 1;
40*54a24607SJunchao Zhang   indices[3] = 1;
41*54a24607SJunchao Zhang   PetscCall(ISCreateGeneral(PETSC_COMM_WORLD,n,indices,PETSC_COPY_VALUES,&is2));
42*54a24607SJunchao Zhang   PetscCall(PetscFree(indices));
43*54a24607SJunchao Zhang   PetscCall(VecCreateMPI(PETSC_COMM_WORLD,n/2,n/2,&vec2));
44*54a24607SJunchao Zhang 
45*54a24607SJunchao Zhang   PetscCall(VecScatterCreate(vec,is,vec2,is2,&scatter));
46*54a24607SJunchao Zhang   PetscCall(ISDestroy(&is));
47*54a24607SJunchao Zhang   PetscCall(ISDestroy(&is2));
48*54a24607SJunchao Zhang 
49*54a24607SJunchao Zhang   PetscCall(VecScatterCopy(scatter,&scatter2));
50*54a24607SJunchao Zhang 
51*54a24607SJunchao Zhang   PetscCall(VecDestroy(&vec));
52*54a24607SJunchao Zhang   PetscCall(VecDestroy(&vec2));
53*54a24607SJunchao Zhang   PetscCall(VecScatterDestroy(&scatter));
54*54a24607SJunchao Zhang   PetscCall(VecScatterDestroy(&scatter2));
55*54a24607SJunchao Zhang   PetscCall(PetscFinalize());
56*54a24607SJunchao Zhang }
57*54a24607SJunchao Zhang 
58*54a24607SJunchao Zhang /*TEST
59*54a24607SJunchao Zhang   test:
60*54a24607SJunchao Zhang     nsize: 1
61*54a24607SJunchao Zhang     output_file: output/empty.out
62*54a24607SJunchao Zhang TEST*/
63