xref: /petsc/src/mat/tests/ex229.c (revision b122ec5aa1bd4469eb4e0673542fb7de3f411254)
1c4762a1bSJed Brown static char help[] = "Test MATMFFD for the rectangular case\n\n";
2c4762a1bSJed Brown 
3c4762a1bSJed Brown #include <petscmat.h>
4c4762a1bSJed Brown 
5c4762a1bSJed Brown static PetscErrorCode myF(void* ctx,Vec x,Vec y)
6c4762a1bSJed Brown {
7c4762a1bSJed Brown   const PetscScalar *ax;
8c4762a1bSJed Brown   PetscScalar       *ay;
9c4762a1bSJed Brown   PetscInt          i,j,m,n;
10c4762a1bSJed Brown 
11c4762a1bSJed Brown   PetscFunctionBegin;
125f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArrayRead(x,&ax));
135f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArray(y,&ay));
145f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetLocalSize(y,&m));
155f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetLocalSize(x,&n));
16c4762a1bSJed Brown   for (i=0;i<m;i++) {
17c4762a1bSJed Brown     PetscScalar xx,yy;
18c4762a1bSJed Brown 
19c4762a1bSJed Brown     yy = 0.0;
20c4762a1bSJed Brown     for (j=0;j<n;j++) {
21c4762a1bSJed Brown       xx = PetscPowScalarInt(ax[j],i+1);
22c4762a1bSJed Brown       yy += xx;
23c4762a1bSJed Brown     }
24c4762a1bSJed Brown     ay[i] = yy;
25c4762a1bSJed Brown   }
265f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArray(y,&ay));
275f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArrayRead(x,&ax));
28c4762a1bSJed Brown   PetscFunctionReturn(0);
29c4762a1bSJed Brown }
30c4762a1bSJed Brown 
31c4762a1bSJed Brown int main(int argc,char **args)
32c4762a1bSJed Brown {
33c4762a1bSJed Brown   Mat            A,B;
34c4762a1bSJed Brown   Vec            base;
35c4762a1bSJed Brown   PetscInt       m = 3 ,n = 2;
36c4762a1bSJed Brown 
37*b122ec5aSJacob Faibussowitsch   CHKERRQ(PetscInitialize(&argc,&args,(char*)0,help));
385f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL));
395f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-m",&m,NULL));
405f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreateMFFD(PETSC_COMM_WORLD,PETSC_DECIDE,PETSC_DECIDE,m,n,&A));
415f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreateVecs(A,&base,NULL));
425f80ce2aSJacob Faibussowitsch   CHKERRQ(VecSet(base,2.0));
435f80ce2aSJacob Faibussowitsch   CHKERRQ(MatMFFDSetFunction(A,myF,NULL));
445f80ce2aSJacob Faibussowitsch   CHKERRQ(MatMFFDSetBase(A,base,NULL));
455f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY));
465f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY));
475f80ce2aSJacob Faibussowitsch   CHKERRQ(MatComputeOperator(A,NULL,&B));
485f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&base));
495f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&A));
505f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&B));
51*b122ec5aSJacob Faibussowitsch   CHKERRQ(PetscFinalize());
52*b122ec5aSJacob Faibussowitsch   return 0;
53c4762a1bSJed Brown }
54c4762a1bSJed Brown 
55c4762a1bSJed Brown /*TEST
56c4762a1bSJed Brown 
57c4762a1bSJed Brown     test:
58c4762a1bSJed Brown       nsize: {{1 2 3 4}}
59c4762a1bSJed Brown 
60c4762a1bSJed Brown TEST*/
61