Actual source code: ex39.c
petsc-3.7.2 2016-06-05
2: static char help[] = "This example is intended for showing how subvectors can\n\
3: share the pointer with the main vector using VecGetArray()\n\
4: and VecPlaceArray() routines so that vector operations done\n\
5: on the subvectors automatically modify the values in the main vector.\n\n";
7: #include <petscvec.h>
9: /* This example shares the array pointers of vectors X,Y,and F with subvectors
10: X1,X2,Y1,Y2,F1,F2 and does vector addition on the subvectors F1 = X1 + Y1, F2 = X2 + Y2 so
11: that F gets updated as a result of sharing the pointers.
12: */
16: int main(int argc,char **argv)
17: {
18: PetscErrorCode ierr;
19: PetscInt N = 10,i;
20: Vec X,Y,F,X1,Y1,X2,Y2,F1,F2;
21: PetscScalar value,zero=0.0;
22: const PetscScalar *x,*y;
23: PetscScalar *f;
25: PetscInitialize(&argc,&argv,(char*)0,help);
27: /* create vectors X,Y and F and set values in it*/
28: VecCreate(PETSC_COMM_SELF,&X);
29: VecSetSizes(X,N,N);
30: VecSetFromOptions(X);
31: VecDuplicate(X,&Y);
32: VecDuplicate(X,&F);
33: PetscObjectSetName((PetscObject)F,"F");
34: for (i=0; i < N; i++) {
35: value = i;
36: VecSetValues(X,1,&i,&value,INSERT_VALUES);
37: value = 100 + i;
38: VecSetValues(Y,1,&i,&value,INSERT_VALUES);
39: }
40: VecSet(F,zero);
42: /* Create subvectors X1,X2,Y1,Y2,F1,F2 */
43: VecCreate(PETSC_COMM_SELF,&X1);
44: VecSetSizes(X1,N/2,N/2);
45: VecSetFromOptions(X1);
46: VecDuplicate(X1,&X2);
47: VecDuplicate(X1,&Y1);
48: VecDuplicate(X1,&Y2);
49: VecDuplicate(X1,&F1);
50: VecDuplicate(X1,&F2);
52: /* Get array pointers for X,Y,F */
53: VecGetArrayRead(X,&x);
54: VecGetArrayRead(Y,&y);
55: VecGetArray(F,&f);
56: /* Share X,Y,F array pointers with subvectors */
57: VecPlaceArray(X1,x);
58: VecPlaceArray(X2,x+N/2);
59: VecPlaceArray(Y1,y);
60: VecPlaceArray(Y2,y+N/2);
61: VecPlaceArray(F1,f);
62: VecPlaceArray(F2,f+N/2);
64: /* Do subvector addition */
65: VecWAXPY(F1,1.0,X1,Y1);
66: VecWAXPY(F2,1.0,X2,Y2);
68: /* Reset subvectors */
69: VecResetArray(X1);
70: VecResetArray(X2);
71: VecResetArray(Y1);
72: VecResetArray(Y2);
73: VecResetArray(F1);
74: VecResetArray(F2);
76: /* Restore X,Y,and F */
77: VecRestoreArrayRead(X,&x);
78: VecRestoreArrayRead(Y,&y);
79: VecRestoreArray(F,&f);
81: PetscPrintf(PETSC_COMM_SELF,"F = X + Y\n");
82: VecView(F,0);
83: /* Destroy vectors */
84: VecDestroy(&X);
85: VecDestroy(&Y);
86: VecDestroy(&F);
87: VecDestroy(&X1);
88: VecDestroy(&Y1);
89: VecDestroy(&F1);
90: VecDestroy(&X2);
91: VecDestroy(&Y2);
92: VecDestroy(&F2);
94: PetscFinalize();
95: return 0;
96: }