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