xref: /petsc/src/mat/tests/ex75.c (revision 327415f76d85372a4417cf1aaa14db707d4d6c04)
1c4762a1bSJed Brown 
2c4762a1bSJed Brown static char help[] = "Tests various routines in MatMPISBAIJ format.\n";
3c4762a1bSJed Brown 
4c4762a1bSJed Brown #include <petscmat.h>
5c4762a1bSJed Brown 
6c4762a1bSJed Brown int main(int argc,char **args)
7c4762a1bSJed Brown {
8c4762a1bSJed Brown   Vec            x,y,u,s1,s2;
9c4762a1bSJed Brown   Mat            A,sA,sB;
10c4762a1bSJed Brown   PetscRandom    rctx;
11c4762a1bSJed Brown   PetscReal      r1,r2,rnorm,tol = PETSC_SQRT_MACHINE_EPSILON;
12c4762a1bSJed Brown   PetscScalar    one=1.0, neg_one=-1.0, value[3], four=4.0,alpha=0.1;
1343359b5eSHong Zhang   PetscInt       n,col[3],n1,block,row,i,j,i2,j2,Ii,J,rstart,rend,bs=1,mbs=16,d_nz=3,o_nz=3,prob=1;
14c4762a1bSJed Brown   PetscMPIInt    size,rank;
15c4762a1bSJed Brown   PetscBool      flg;
16c4762a1bSJed Brown 
17*327415f7SBarry Smith   PetscFunctionBeginUser;
189566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc,&args,(char*)0,help));
199566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL,NULL,"-mbs",&mbs,NULL));
209566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL,NULL,"-bs",&bs,NULL));
219566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL,NULL,"-prob",&prob,NULL));
22c4762a1bSJed Brown 
239566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&rank));
249566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size));
25c4762a1bSJed Brown 
2643359b5eSHong Zhang   /* Create a BAIJ matrix A */
27c4762a1bSJed Brown   n = mbs*bs;
289566063dSJacob Faibussowitsch   PetscCall(MatCreate(PETSC_COMM_WORLD,&A));
299566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,n,n));
309566063dSJacob Faibussowitsch   PetscCall(MatSetType(A,MATBAIJ));
319566063dSJacob Faibussowitsch   PetscCall(MatSetFromOptions(A));
329566063dSJacob Faibussowitsch   PetscCall(MatMPIBAIJSetPreallocation(A,bs,d_nz,NULL,o_nz,NULL));
339566063dSJacob Faibussowitsch   PetscCall(MatSeqBAIJSetPreallocation(A,bs,d_nz,NULL));
349566063dSJacob Faibussowitsch   PetscCall(MatSetOption(A,MAT_NEW_NONZERO_ALLOCATION_ERR,PETSC_FALSE));
35c4762a1bSJed Brown 
36c4762a1bSJed Brown   if (bs == 1) {
37c4762a1bSJed Brown     if (prob == 1) { /* tridiagonal matrix */
3843359b5eSHong Zhang       value[0] = -1.0; value[1] = 0.0; value[2] = -1.0;
39c4762a1bSJed Brown       for (i=1; i<n-1; i++) {
40c4762a1bSJed Brown         col[0] = i-1; col[1] = i; col[2] = i+1;
419566063dSJacob Faibussowitsch         PetscCall(MatSetValues(A,1,&i,3,col,value,INSERT_VALUES));
42c4762a1bSJed Brown       }
43c4762a1bSJed Brown       i       = n - 1; col[0]=0; col[1] = n - 2; col[2] = n - 1;
4443359b5eSHong Zhang       value[0]= 0.1; value[1]=-1.0; value[2]=0.0;
459566063dSJacob Faibussowitsch       PetscCall(MatSetValues(A,1,&i,3,col,value,INSERT_VALUES));
46c4762a1bSJed Brown 
47c4762a1bSJed Brown       i        = 0; col[0] = 0; col[1] = 1; col[2]=n-1;
4843359b5eSHong Zhang       value[0] = 0.0; value[1] = -1.0; value[2]=0.1;
499566063dSJacob Faibussowitsch       PetscCall(MatSetValues(A,1,&i,3,col,value,INSERT_VALUES));
50c4762a1bSJed Brown     } else if (prob ==2) { /* matrix for the five point stencil */
51c4762a1bSJed Brown       n1 = (int) PetscSqrtReal((PetscReal)n);
52c4762a1bSJed Brown       for (i=0; i<n1; i++) {
53c4762a1bSJed Brown         for (j=0; j<n1; j++) {
54c4762a1bSJed Brown           Ii = j + n1*i;
559566063dSJacob Faibussowitsch           if (i>0)    {J = Ii - n1; PetscCall(MatSetValues(A,1,&Ii,1,&J,&neg_one,INSERT_VALUES));}
569566063dSJacob Faibussowitsch           if (i<n1-1) {J = Ii + n1; PetscCall(MatSetValues(A,1,&Ii,1,&J,&neg_one,INSERT_VALUES));}
579566063dSJacob Faibussowitsch           if (j>0)    {J = Ii - 1;  PetscCall(MatSetValues(A,1,&Ii,1,&J,&neg_one,INSERT_VALUES));}
589566063dSJacob Faibussowitsch           if (j<n1-1) {J = Ii + 1;  PetscCall(MatSetValues(A,1,&Ii,1,&J,&neg_one,INSERT_VALUES));}
599566063dSJacob Faibussowitsch           PetscCall(MatSetValues(A,1,&Ii,1,&Ii,&four,INSERT_VALUES));
60c4762a1bSJed Brown         }
61c4762a1bSJed Brown       }
62c4762a1bSJed Brown     }
63c4762a1bSJed Brown     /* end of if (bs == 1) */
64c4762a1bSJed Brown   } else {  /* bs > 1 */
65c4762a1bSJed Brown     for (block=0; block<n/bs; block++) {
66c4762a1bSJed Brown       /* diagonal blocks */
67c4762a1bSJed Brown       value[0] = -1.0; value[1] = 4.0; value[2] = -1.0;
68c4762a1bSJed Brown       for (i=1+block*bs; i<bs-1+block*bs; i++) {
69c4762a1bSJed Brown         col[0] = i-1; col[1] = i; col[2] = i+1;
709566063dSJacob Faibussowitsch         PetscCall(MatSetValues(A,1,&i,3,col,value,INSERT_VALUES));
71c4762a1bSJed Brown       }
72c4762a1bSJed Brown       i       = bs - 1+block*bs; col[0] = bs - 2+block*bs; col[1] = bs - 1+block*bs;
73c4762a1bSJed Brown       value[0]=-1.0; value[1]=4.0;
749566063dSJacob Faibussowitsch       PetscCall(MatSetValues(A,1,&i,2,col,value,INSERT_VALUES));
75c4762a1bSJed Brown 
76c4762a1bSJed Brown       i       = 0+block*bs; col[0] = 0+block*bs; col[1] = 1+block*bs;
77c4762a1bSJed Brown       value[0]=4.0; value[1] = -1.0;
789566063dSJacob Faibussowitsch       PetscCall(MatSetValues(A,1,&i,2,col,value,INSERT_VALUES));
79c4762a1bSJed Brown     }
80c4762a1bSJed Brown     /* off-diagonal blocks */
81c4762a1bSJed Brown     value[0]=-1.0;
82c4762a1bSJed Brown     for (i=0; i<(n/bs-1)*bs; i++) {
83c4762a1bSJed Brown       col[0]=i+bs;
849566063dSJacob Faibussowitsch       PetscCall(MatSetValues(A,1,&i,1,col,value,INSERT_VALUES));
85c4762a1bSJed Brown       col[0]=i; row=i+bs;
869566063dSJacob Faibussowitsch       PetscCall(MatSetValues(A,1,&row,1,col,value,INSERT_VALUES));
87c4762a1bSJed Brown     }
88c4762a1bSJed Brown   }
899566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY));
909566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY));
919566063dSJacob Faibussowitsch   PetscCall(MatSetOption(A,MAT_SYMMETRIC,PETSC_TRUE));
9243359b5eSHong Zhang 
9343359b5eSHong Zhang   /* Get SBAIJ matrix sA from A */
949566063dSJacob Faibussowitsch   PetscCall(MatConvert(A,MATSBAIJ,MAT_INITIAL_MATRIX,&sA));
95c4762a1bSJed Brown 
96c4762a1bSJed Brown   /* Test MatGetSize(), MatGetLocalSize() */
979566063dSJacob Faibussowitsch   PetscCall(MatGetSize(sA, &i,&j));
989566063dSJacob Faibussowitsch   PetscCall(MatGetSize(A, &i2,&j2));
99c4762a1bSJed Brown   i   -= i2; j -= j2;
100c4762a1bSJed Brown   if (i || j) {
1019566063dSJacob Faibussowitsch     PetscCall(PetscSynchronizedPrintf(PETSC_COMM_WORLD,"[%d], Error: MatGetSize()\n",rank));
1029566063dSJacob Faibussowitsch     PetscCall(PetscSynchronizedFlush(PETSC_COMM_WORLD,PETSC_STDOUT));
103c4762a1bSJed Brown   }
104c4762a1bSJed Brown 
1059566063dSJacob Faibussowitsch   PetscCall(MatGetLocalSize(sA, &i,&j));
1069566063dSJacob Faibussowitsch   PetscCall(MatGetLocalSize(A, &i2,&j2));
107c4762a1bSJed Brown   i2  -= i; j2 -= j;
108c4762a1bSJed Brown   if (i2 || j2) {
1099566063dSJacob Faibussowitsch     PetscCall(PetscSynchronizedPrintf(PETSC_COMM_WORLD,"[%d], Error: MatGetLocalSize()\n",rank));
1109566063dSJacob Faibussowitsch     PetscCall(PetscSynchronizedFlush(PETSC_COMM_WORLD,PETSC_STDOUT));
111c4762a1bSJed Brown   }
112c4762a1bSJed Brown 
113c4762a1bSJed Brown   /* vectors */
114c4762a1bSJed Brown   /*--------------------*/
115c4762a1bSJed Brown   /* i is obtained from MatGetLocalSize() */
1169566063dSJacob Faibussowitsch   PetscCall(VecCreate(PETSC_COMM_WORLD,&x));
1179566063dSJacob Faibussowitsch   PetscCall(VecSetSizes(x,i,PETSC_DECIDE));
1189566063dSJacob Faibussowitsch   PetscCall(VecSetFromOptions(x));
1199566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(x,&y));
1209566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(x,&u));
1219566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(x,&s1));
1229566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(x,&s2));
123c4762a1bSJed Brown 
1249566063dSJacob Faibussowitsch   PetscCall(PetscRandomCreate(PETSC_COMM_WORLD,&rctx));
1259566063dSJacob Faibussowitsch   PetscCall(PetscRandomSetFromOptions(rctx));
1269566063dSJacob Faibussowitsch   PetscCall(VecSetRandom(x,rctx));
1279566063dSJacob Faibussowitsch   PetscCall(VecSet(u,one));
128c4762a1bSJed Brown 
129c4762a1bSJed Brown   /* Test MatNorm() */
1309566063dSJacob Faibussowitsch   PetscCall(MatNorm(A,NORM_FROBENIUS,&r1));
1319566063dSJacob Faibussowitsch   PetscCall(MatNorm(sA,NORM_FROBENIUS,&r2));
132c4762a1bSJed Brown   rnorm = PetscAbsReal(r1-r2)/r2;
133dd400576SPatrick Sanan   if (rnorm > tol && rank == 0) {
1349566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(PETSC_COMM_SELF,"Error: MatNorm_FROBENIUS(), Anorm=%16.14e, sAnorm=%16.14e bs=%" PetscInt_FMT "\n",(double)r1,(double)r2,bs));
135c4762a1bSJed Brown   }
1369566063dSJacob Faibussowitsch   PetscCall(MatNorm(A,NORM_INFINITY,&r1));
1379566063dSJacob Faibussowitsch   PetscCall(MatNorm(sA,NORM_INFINITY,&r2));
138c4762a1bSJed Brown   rnorm = PetscAbsReal(r1-r2)/r2;
139dd400576SPatrick Sanan   if (rnorm > tol && rank == 0) {
1409566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Error: MatNorm_INFINITY(), Anorm=%16.14e, sAnorm=%16.14e bs=%" PetscInt_FMT "\n",(double)r1,(double)r2,bs));
141c4762a1bSJed Brown   }
1429566063dSJacob Faibussowitsch   PetscCall(MatNorm(A,NORM_1,&r1));
1439566063dSJacob Faibussowitsch   PetscCall(MatNorm(sA,NORM_1,&r2));
144c4762a1bSJed Brown   rnorm = PetscAbsReal(r1-r2)/r2;
145dd400576SPatrick Sanan   if (rnorm > tol && rank == 0) {
1469566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Error: MatNorm_1(), Anorm=%16.14e, sAnorm=%16.14e bs=%" PetscInt_FMT "\n",(double)r1,(double)r2,bs));
147c4762a1bSJed Brown   }
148c4762a1bSJed Brown 
149c4762a1bSJed Brown   /* Test MatGetOwnershipRange() */
1509566063dSJacob Faibussowitsch   PetscCall(MatGetOwnershipRange(sA,&rstart,&rend));
1519566063dSJacob Faibussowitsch   PetscCall(MatGetOwnershipRange(A,&i2,&j2));
152c4762a1bSJed Brown   i2  -= rstart; j2 -= rend;
153c4762a1bSJed Brown   if (i2 || j2) {
1549566063dSJacob Faibussowitsch     PetscCall(PetscSynchronizedPrintf(PETSC_COMM_WORLD,"[%d], Error: MaGetOwnershipRange()\n",rank));
1559566063dSJacob Faibussowitsch     PetscCall(PetscSynchronizedFlush(PETSC_COMM_WORLD,PETSC_STDOUT));
156c4762a1bSJed Brown   }
157c4762a1bSJed Brown 
158c4762a1bSJed Brown   /* Test MatDiagonalScale() */
1599566063dSJacob Faibussowitsch   PetscCall(MatDiagonalScale(A,x,x));
1609566063dSJacob Faibussowitsch   PetscCall(MatDiagonalScale(sA,x,x));
1619566063dSJacob Faibussowitsch   PetscCall(MatMultEqual(A,sA,10,&flg));
16228b400f6SJacob Faibussowitsch   PetscCheck(flg,PETSC_COMM_SELF,PETSC_ERR_ARG_NOTSAMETYPE,"Error in MatDiagonalScale");
163c4762a1bSJed Brown 
164c4762a1bSJed Brown   /* Test MatGetDiagonal(), MatScale() */
1659566063dSJacob Faibussowitsch   PetscCall(MatGetDiagonal(A,s1));
1669566063dSJacob Faibussowitsch   PetscCall(MatGetDiagonal(sA,s2));
1679566063dSJacob Faibussowitsch   PetscCall(VecNorm(s1,NORM_1,&r1));
1689566063dSJacob Faibussowitsch   PetscCall(VecNorm(s2,NORM_1,&r2));
169c4762a1bSJed Brown   r1  -= r2;
170c4762a1bSJed Brown   if (r1<-tol || r1>tol) {
1719566063dSJacob Faibussowitsch     PetscCall(PetscSynchronizedPrintf(PETSC_COMM_WORLD,"[%d], Error: MatDiagonalScale() or MatGetDiagonal(), r1=%g \n",rank,(double)r1));
1729566063dSJacob Faibussowitsch     PetscCall(PetscSynchronizedFlush(PETSC_COMM_WORLD,PETSC_STDOUT));
173c4762a1bSJed Brown   }
174c4762a1bSJed Brown 
1759566063dSJacob Faibussowitsch   PetscCall(MatScale(A,alpha));
1769566063dSJacob Faibussowitsch   PetscCall(MatScale(sA,alpha));
177c4762a1bSJed Brown 
178c4762a1bSJed Brown   /* Test MatGetRowMaxAbs() */
1799566063dSJacob Faibussowitsch   PetscCall(MatGetRowMaxAbs(A,s1,NULL));
1809566063dSJacob Faibussowitsch   PetscCall(MatGetRowMaxAbs(sA,s2,NULL));
181c4762a1bSJed Brown 
1829566063dSJacob Faibussowitsch   PetscCall(VecNorm(s1,NORM_1,&r1));
1839566063dSJacob Faibussowitsch   PetscCall(VecNorm(s2,NORM_1,&r2));
184c4762a1bSJed Brown   r1  -= r2;
185c4762a1bSJed Brown   if (r1<-tol || r1>tol) {
1869566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(PETSC_COMM_SELF,"Error: MatGetRowMaxAbs() \n"));
187c4762a1bSJed Brown   }
188c4762a1bSJed Brown 
189c4762a1bSJed Brown   /* Test MatMult(), MatMultAdd() */
1909566063dSJacob Faibussowitsch   PetscCall(MatMultEqual(A,sA,10,&flg));
191c4762a1bSJed Brown   if (!flg) {
1929566063dSJacob Faibussowitsch     PetscCall(PetscSynchronizedPrintf(PETSC_COMM_WORLD,"[%d], Error: MatMult() or MatScale()\n",rank));
1939566063dSJacob Faibussowitsch     PetscCall(PetscSynchronizedFlush(PETSC_COMM_WORLD,PETSC_STDOUT));
194c4762a1bSJed Brown   }
195c4762a1bSJed Brown 
1969566063dSJacob Faibussowitsch   PetscCall(MatMultAddEqual(A,sA,10,&flg));
197c4762a1bSJed Brown   if (!flg) {
1989566063dSJacob Faibussowitsch     PetscCall(PetscSynchronizedPrintf(PETSC_COMM_WORLD,"[%d], Error: MatMultAdd()\n",rank));
1999566063dSJacob Faibussowitsch     PetscCall(PetscSynchronizedFlush(PETSC_COMM_WORLD,PETSC_STDOUT));
200c4762a1bSJed Brown   }
201c4762a1bSJed Brown 
202c4762a1bSJed Brown   /* Test MatMultTranspose(), MatMultTransposeAdd() */
203c4762a1bSJed Brown   for (i=0; i<10; i++) {
2049566063dSJacob Faibussowitsch     PetscCall(VecSetRandom(x,rctx));
2059566063dSJacob Faibussowitsch     PetscCall(MatMultTranspose(A,x,s1));
2069566063dSJacob Faibussowitsch     PetscCall(MatMultTranspose(sA,x,s2));
2079566063dSJacob Faibussowitsch     PetscCall(VecNorm(s1,NORM_1,&r1));
2089566063dSJacob Faibussowitsch     PetscCall(VecNorm(s2,NORM_1,&r2));
209c4762a1bSJed Brown     r1  -= r2;
210c4762a1bSJed Brown     if (r1<-tol || r1>tol) {
2119566063dSJacob Faibussowitsch       PetscCall(PetscSynchronizedPrintf(PETSC_COMM_WORLD,"[%d], Error: MatMult() or MatScale(), err=%g\n",rank,(double)r1));
2129566063dSJacob Faibussowitsch       PetscCall(PetscSynchronizedFlush(PETSC_COMM_WORLD,PETSC_STDOUT));
213c4762a1bSJed Brown     }
214c4762a1bSJed Brown   }
215c4762a1bSJed Brown   for (i=0; i<10; i++) {
2169566063dSJacob Faibussowitsch     PetscCall(VecSetRandom(x,rctx));
2179566063dSJacob Faibussowitsch     PetscCall(VecSetRandom(y,rctx));
2189566063dSJacob Faibussowitsch     PetscCall(MatMultTransposeAdd(A,x,y,s1));
2199566063dSJacob Faibussowitsch     PetscCall(MatMultTransposeAdd(sA,x,y,s2));
2209566063dSJacob Faibussowitsch     PetscCall(VecNorm(s1,NORM_1,&r1));
2219566063dSJacob Faibussowitsch     PetscCall(VecNorm(s2,NORM_1,&r2));
222c4762a1bSJed Brown     r1  -= r2;
223c4762a1bSJed Brown     if (r1<-tol || r1>tol) {
2249566063dSJacob Faibussowitsch       PetscCall(PetscSynchronizedPrintf(PETSC_COMM_WORLD,"[%d], Error: MatMultAdd(), err=%g \n",rank,(double)r1));
2259566063dSJacob Faibussowitsch       PetscCall(PetscSynchronizedFlush(PETSC_COMM_WORLD,PETSC_STDOUT));
226c4762a1bSJed Brown     }
227c4762a1bSJed Brown   }
228c4762a1bSJed Brown 
229c4762a1bSJed Brown   /* Test MatDuplicate() */
2309566063dSJacob Faibussowitsch   PetscCall(MatDuplicate(sA,MAT_COPY_VALUES,&sB));
2319566063dSJacob Faibussowitsch   PetscCall(MatEqual(sA,sB,&flg));
232c4762a1bSJed Brown   if (!flg) {
2339566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(PETSC_COMM_WORLD," Error in MatDuplicate(), sA != sB \n"));
234c4762a1bSJed Brown   }
2359566063dSJacob Faibussowitsch   PetscCall(MatMultEqual(sA,sB,5,&flg));
236c4762a1bSJed Brown   if (!flg) {
2379566063dSJacob Faibussowitsch     PetscCall(PetscSynchronizedPrintf(PETSC_COMM_WORLD,"[%d], Error: MatDuplicate() or MatMult()\n",rank));
2389566063dSJacob Faibussowitsch     PetscCall(PetscSynchronizedFlush(PETSC_COMM_WORLD,PETSC_STDOUT));
239c4762a1bSJed Brown   }
2409566063dSJacob Faibussowitsch   PetscCall(MatMultAddEqual(sA,sB,5,&flg));
241c4762a1bSJed Brown   if (!flg) {
2429566063dSJacob Faibussowitsch     PetscCall(PetscSynchronizedPrintf(PETSC_COMM_WORLD,"[%d], Error: MatDuplicate() or MatMultAdd(()\n",rank));
2439566063dSJacob Faibussowitsch     PetscCall(PetscSynchronizedFlush(PETSC_COMM_WORLD,PETSC_STDOUT));
244c4762a1bSJed Brown   }
2459566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&sB));
2469566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&u));
2479566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&x));
2489566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&y));
2499566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&s1));
2509566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&s2));
2519566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&sA));
2529566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&A));
2539566063dSJacob Faibussowitsch   PetscCall(PetscRandomDestroy(&rctx));
2549566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
255b122ec5aSJacob Faibussowitsch   return 0;
256c4762a1bSJed Brown }
257c4762a1bSJed Brown 
258c4762a1bSJed Brown /*TEST
259c4762a1bSJed Brown 
260c4762a1bSJed Brown    test:
261c4762a1bSJed Brown       nsize: {{1 3}}
26243359b5eSHong Zhang       args: -bs {{1 2 3  5  7 8}} -mat_ignore_lower_triangular -prob {{1 2}}
263c4762a1bSJed Brown 
264c4762a1bSJed Brown TEST*/
265