Actual source code: ex28.c
2: static char help[] = "Tests repeated VecDotBegin()/VecDotEnd().\n\n";
4: #include petscvec.h
8: int main(int argc,char **argv)
9: {
11: PetscInt n = 25,i,row0 = 0;
12: PetscScalar one = 1.0,two = 2.0,result1,result2,results[40],value,ten = 10.0;
13: PetscScalar result1a,result2a;
14: PetscReal result3,result4,result[2],result3a,result4a,resulta[2];
15: Vec x,y,vecs[40];
17: PetscInitialize(&argc,&argv,(char*)0,help);
19: /* create vector */
20: VecCreate(PETSC_COMM_WORLD,&x);
21: VecSetSizes(x,n,PETSC_DECIDE);
22: VecSetFromOptions(x);
23: VecDuplicate(x,&y);
25: VecSet(x,one);
26: VecSet(y,two);
28: /*
29: Test mixing dot products and norms that require sums
30: */
31: result1 = result2 = 0.0;
32: result3 = result4 = 0.0;
33: VecDotBegin(x,y,&result1);
34: VecDotBegin(y,x,&result2);
35: VecNormBegin(y,NORM_2,&result3);
36: VecNormBegin(x,NORM_1,&result4);
37: VecDotEnd(x,y,&result1);
38: VecDotEnd(y,x,&result2);
39: VecNormEnd(y,NORM_2,&result3);
40: VecNormEnd(x,NORM_1,&result4);
41:
42: VecDot(x,y,&result1a);
43: VecDot(y,x,&result2a);
44: VecNorm(y,NORM_2,&result3a);
45: VecNorm(x,NORM_1,&result4a);
46:
47: if (result1 != result1a || result2 != result2a) {
48: PetscPrintf(PETSC_COMM_WORLD,"Error dot: result1 %G result2 %G\n",PetscRealPart(result1),PetscRealPart(result2));
49: }
50: if (result3 != result3a || result4 != result4a) {
51: PetscPrintf(PETSC_COMM_WORLD,"Error 1,2 norms: result3 %G result4 %G\n",result3,result4);
52: }
54: /*
55: Test norms that only require abs
56: */
57: result1 = result2 = 0.0;
58: result3 = result4 = 0.0;
59: VecNormBegin(y,NORM_MAX,&result3);
60: VecNormBegin(x,NORM_MAX,&result4);
61: VecNormEnd(y,NORM_MAX,&result3);
62: VecNormEnd(x,NORM_MAX,&result4);
64: VecNorm(x,NORM_MAX,&result4a);
65: VecNorm(y,NORM_MAX,&result3a);
66: if (result3 != result3a || result4 != result4a) {
67: PetscPrintf(PETSC_COMM_WORLD,"Error max norm: result3 %G result4 %G\n",result3,result4);
68: }
70: /*
71: Tests dot, max, 1, norm
72: */
73: result1 = result2 = 0.0;
74: result3 = result4 = 0.0;
75: VecSetValues(x,1,&row0,&ten,INSERT_VALUES);
76: VecAssemblyBegin(x);
77: VecAssemblyEnd(x);
79: VecDotBegin(x,y,&result1);
80: VecDotBegin(y,x,&result2);
81: VecNormBegin(x,NORM_MAX,&result3);
82: VecNormBegin(x,NORM_1,&result4);
83: VecDotEnd(x,y,&result1);
84: VecDotEnd(y,x,&result2);
85: VecNormEnd(x,NORM_MAX,&result3);
86: VecNormEnd(x,NORM_1,&result4);
88: VecDot(x,y,&result1a);
89: VecDot(y,x,&result2a);
90: VecNorm(x,NORM_MAX,&result3a);
91: VecNorm(x,NORM_1,&result4a);
92: if (result1 != result1a || result2 != result2a) {
93: PetscPrintf(PETSC_COMM_WORLD,"Error dot: result1 %G result2 %G\n",PetscRealPart(result1),PetscRealPart(result2));
94: }
95: if (result3 != result3a || result4 != result4a) {
96: PetscPrintf(PETSC_COMM_WORLD,"Error max 1 norms: result3 %G result4 %G\n",result3,result4);
97: }
99: /*
100: tests 1_and_2 norm
101: */
102: VecNormBegin(x,NORM_MAX,&result3);
103: VecNormBegin(x,NORM_1_AND_2,result);
104: VecNormBegin(y,NORM_MAX,&result4);
105: VecNormEnd(x,NORM_MAX,&result3);
106: VecNormEnd(x,NORM_1_AND_2,result);
107: VecNormEnd(y,NORM_MAX,&result4);
109: VecNorm(x,NORM_MAX,&result3a);
110: VecNorm(x,NORM_1_AND_2,resulta);
111: VecNorm(y,NORM_MAX,&result4a);
112: if (result3 != result3a || result4 != result4a) {
113: PetscPrintf(PETSC_COMM_WORLD,"Error max: result1 %G result2 %G\n",result3,result4);
114: }
115: if (PetscAbsReal(result[0]-resulta[0]) > .01 || PetscAbsReal(result[1]-resulta[1]) > .01) {
116: PetscPrintf(PETSC_COMM_WORLD,"Error 1 and 2 norms: result[0] %G result[1] %G\n",result[0],result[1]);
117: }
119: VecDestroy(x);
120: VecDestroy(y);
122: /*
123: Tests computing a large number of operations that require
124: allocating a larger data structure internally
125: */
126: for (i=0; i<40; i++) {
127: VecCreate(PETSC_COMM_WORLD,vecs+i);
128: VecSetSizes(vecs[i],PETSC_DECIDE,n);
129: VecSetFromOptions(vecs[i]);
130: value = (PetscReal)i;
131: VecSet(vecs[i],value);
132: }
133: for (i=0; i<39; i++) {
134: VecDotBegin(vecs[i],vecs[i+1],results+i);
135: }
136: for (i=0; i<39; i++) {
137: VecDotEnd(vecs[i],vecs[i+1],results+i);
138: if (results[i] != 25.0*i*(i+1)) {
139: PetscPrintf(PETSC_COMM_WORLD,"i %D expected %G got %G\n",i,25.0*i*(i+1),PetscRealPart(results[i]));
140: }
141: }
142: for (i=0; i<40; i++) {
143: VecDestroy(vecs[i]);
144: }
146: PetscFinalize();
147: return 0;
148: }
149: