xref: /petsc/src/mat/tests/ex229.c (revision ccb4e88a40f0b86eaeca07ff64c64e4de2fae686)
1 static char help[] = "Test MATMFFD for the rectangular case\n\n";
2 
3 #include <petscmat.h>
4 
5 static PetscErrorCode myF(void* ctx,Vec x,Vec y)
6 {
7   const PetscScalar *ax;
8   PetscScalar       *ay;
9   PetscInt          i,j,m,n;
10   PetscErrorCode    ierr;
11 
12   PetscFunctionBegin;
13   ierr = VecGetArrayRead(x,&ax);CHKERRQ(ierr);
14   ierr = VecGetArray(y,&ay);CHKERRQ(ierr);
15   ierr = VecGetLocalSize(y,&m);CHKERRQ(ierr);
16   ierr = VecGetLocalSize(x,&n);CHKERRQ(ierr);
17   for (i=0;i<m;i++) {
18     PetscScalar xx,yy;
19 
20     yy = 0.0;
21     for (j=0;j<n;j++) {
22       xx = PetscPowScalarInt(ax[j],i+1);
23       yy += xx;
24     }
25     ay[i] = yy;
26   }
27   ierr = VecRestoreArray(y,&ay);CHKERRQ(ierr);
28   ierr = VecRestoreArrayRead(x,&ax);CHKERRQ(ierr);
29   PetscFunctionReturn(0);
30 }
31 
32 int main(int argc,char **args)
33 {
34   Mat            A,B;
35   Vec            base;
36   PetscInt       m = 3 ,n = 2;
37   PetscErrorCode ierr;
38 
39   ierr = PetscInitialize(&argc,&args,(char*)0,help);if (ierr) return ierr;
40   ierr = PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL);CHKERRQ(ierr);
41   ierr = PetscOptionsGetInt(NULL,NULL,"-m",&m,NULL);CHKERRQ(ierr);
42   ierr = MatCreateMFFD(PETSC_COMM_WORLD,PETSC_DECIDE,PETSC_DECIDE,m,n,&A);CHKERRQ(ierr);
43   ierr = MatCreateVecs(A,&base,NULL);CHKERRQ(ierr);
44   ierr = VecSet(base,2.0);CHKERRQ(ierr);
45   ierr = MatMFFDSetFunction(A,myF,NULL);CHKERRQ(ierr);
46   ierr = MatMFFDSetBase(A,base,NULL);CHKERRQ(ierr);
47   ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
48   ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
49   ierr = MatComputeOperator(A,NULL,&B);CHKERRQ(ierr);
50   ierr = VecDestroy(&base);CHKERRQ(ierr);
51   ierr = MatDestroy(&A);CHKERRQ(ierr);
52   ierr = MatDestroy(&B);CHKERRQ(ierr);
53   ierr = PetscFinalize();
54   return ierr;
55 }
56 
57 /*TEST
58 
59     test:
60       nsize: {{1 2 3 4}}
61 
62 TEST*/
63