1c4762a1bSJed Brown static char help[] = "Tests MAIJ matrix for large DOF\n\n"; 2c4762a1bSJed Brown 3c4762a1bSJed Brown #include <petscdm.h> 4c4762a1bSJed Brown #include <petscdmda.h> 5c4762a1bSJed Brown 6c4762a1bSJed Brown int main(int argc,char *argv[]) 7c4762a1bSJed Brown { 8c4762a1bSJed Brown Mat M; 9c4762a1bSJed Brown Vec x,y; 10c4762a1bSJed Brown PetscErrorCode ierr; 11c4762a1bSJed Brown DM da,daf; 12c4762a1bSJed Brown 13c4762a1bSJed Brown ierr = PetscInitialize(&argc,&argv,0,help);if (ierr) return ierr; 14*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMDACreate2d(PETSC_COMM_WORLD, DM_BOUNDARY_NONE, DM_BOUNDARY_NONE,DMDA_STENCIL_STAR,4,5,PETSC_DECIDE,PETSC_DECIDE,41,1,0,0,&da)); 15*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetFromOptions(da)); 16*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetUp(da)); 17*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMRefine(da,PETSC_COMM_WORLD,&daf)); 18*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMCreateInterpolation(da,daf,&M,NULL)); 19*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMCreateGlobalVector(da,&x)); 20*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMCreateGlobalVector(daf,&y)); 21c4762a1bSJed Brown 22*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatMult(M,x,y)); 23*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatMultTranspose(M,y,x)); 24*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMDestroy(&da)); 25*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMDestroy(&daf)); 26*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&x)); 27*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&y)); 28*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&M)); 29c4762a1bSJed Brown ierr = PetscFinalize(); 30c4762a1bSJed Brown return ierr; 31c4762a1bSJed Brown } 32c4762a1bSJed Brown 33c4762a1bSJed Brown /*TEST 34c4762a1bSJed Brown 35c4762a1bSJed Brown test: 36c4762a1bSJed Brown nsize: 2 37c4762a1bSJed Brown 38c4762a1bSJed Brown TEST*/ 39