Actual source code: ex39.c
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: */
14: int main(int argc,char **argv)
15: {
16: PetscInt N = 10,i;
17: Vec X,Y,F,X1,Y1,X2,Y2,F1,F2;
18: PetscScalar value,zero=0.0;
19: const PetscScalar *x,*y;
20: PetscScalar *f;
22: PetscInitialize(&argc,&argv,(char*)0,help);
24: /* create vectors X,Y and F and set values in it*/
25: VecCreate(PETSC_COMM_SELF,&X);
26: VecSetSizes(X,N,N);
27: VecSetFromOptions(X);
28: VecDuplicate(X,&Y);
29: VecDuplicate(X,&F);
30: PetscObjectSetName((PetscObject)F,"F");
31: for (i=0; i < N; i++) {
32: value = i;
33: VecSetValues(X,1,&i,&value,INSERT_VALUES);
34: value = 100 + i;
35: VecSetValues(Y,1,&i,&value,INSERT_VALUES);
36: }
37: VecSet(F,zero);
39: /* Create subvectors X1,X2,Y1,Y2,F1,F2 */
40: VecCreate(PETSC_COMM_SELF,&X1);
41: VecSetSizes(X1,N/2,N/2);
42: VecSetFromOptions(X1);
43: VecDuplicate(X1,&X2);
44: VecDuplicate(X1,&Y1);
45: VecDuplicate(X1,&Y2);
46: VecDuplicate(X1,&F1);
47: VecDuplicate(X1,&F2);
49: /* Get array pointers for X,Y,F */
50: VecGetArrayRead(X,&x);
51: VecGetArrayRead(Y,&y);
52: VecGetArray(F,&f);
53: /* Share X,Y,F array pointers with subvectors */
54: VecPlaceArray(X1,x);
55: VecPlaceArray(X2,x+N/2);
56: VecPlaceArray(Y1,y);
57: VecPlaceArray(Y2,y+N/2);
58: VecPlaceArray(F1,f);
59: VecPlaceArray(F2,f+N/2);
61: /* Do subvector addition */
62: VecWAXPY(F1,1.0,X1,Y1);
63: VecWAXPY(F2,1.0,X2,Y2);
65: /* Reset subvectors */
66: VecResetArray(X1);
67: VecResetArray(X2);
68: VecResetArray(Y1);
69: VecResetArray(Y2);
70: VecResetArray(F1);
71: VecResetArray(F2);
73: /* Restore X,Y,and F */
74: VecRestoreArrayRead(X,&x);
75: VecRestoreArrayRead(Y,&y);
76: VecRestoreArray(F,&f);
78: PetscPrintf(PETSC_COMM_SELF,"F = X + Y\n");
79: VecView(F,0);
80: /* Destroy vectors */
81: VecDestroy(&X);
82: VecDestroy(&Y);
83: VecDestroy(&F);
84: VecDestroy(&X1);
85: VecDestroy(&Y1);
86: VecDestroy(&F1);
87: VecDestroy(&X2);
88: VecDestroy(&Y2);
89: VecDestroy(&F2);
91: PetscFinalize();
92: return 0;
93: }