1c4762a1bSJed Brown static char help[] = "Tests the use of MatZeroRowsColumns() for parallel matrices.\n\ 2c4762a1bSJed Brown Contributed-by: Stephan Kramer <s.kramer@imperial.ac.uk>\n\n"; 3c4762a1bSJed Brown 4c4762a1bSJed Brown #include <petscmat.h> 5c4762a1bSJed Brown 6c4762a1bSJed Brown int main(int argc,char **args) 7c4762a1bSJed Brown { 8c4762a1bSJed Brown Mat A; 9c4762a1bSJed Brown Vec x, rhs, y; 10c4762a1bSJed Brown PetscInt i,j,k,b,m = 3,n,nlocal=2,bs=1,Ii,J; 11c4762a1bSJed Brown PetscInt *boundary_nodes, nboundary_nodes, *boundary_indices; 12c4762a1bSJed Brown PetscMPIInt rank,size; 13c4762a1bSJed Brown PetscErrorCode ierr; 14c4762a1bSJed Brown PetscScalar v,v0,v1,v2,a0=0.1,a,rhsval, *boundary_values,diag = 1.0; 15c4762a1bSJed Brown PetscReal norm; 16c4762a1bSJed Brown char convname[64]; 17c4762a1bSJed Brown PetscBool upwind = PETSC_FALSE, nonlocalBC = PETSC_FALSE, zerorhs = PETSC_TRUE, convert = PETSC_FALSE; 18c4762a1bSJed Brown 19c4762a1bSJed Brown ierr = PetscInitialize(&argc,&args,(char*)0,help);if (ierr) return ierr; 20*5f80ce2aSJacob Faibussowitsch CHKERRMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&rank)); 21*5f80ce2aSJacob Faibussowitsch CHKERRMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size)); 22c4762a1bSJed Brown n = nlocal*size; 23c4762a1bSJed Brown 24*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetInt(NULL,NULL, "-bs", &bs, NULL)); 25*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetBool(NULL,NULL, "-nonlocal_bc", &nonlocalBC, NULL)); 26*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetScalar(NULL,NULL, "-diag", &diag, NULL)); 27*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetString(NULL,NULL,"-convname",convname,sizeof(convname),&convert)); 28*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetBool(NULL,NULL, "-zerorhs", &zerorhs, NULL)); 29c4762a1bSJed Brown 30*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreate(PETSC_COMM_WORLD,&A)); 31*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,m*n*bs,m*n*bs)); 32*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetFromOptions(A)); 33*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetUp(A)); 34c4762a1bSJed Brown 35*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreateVecs(A, NULL, &rhs)); 36*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecSetFromOptions(rhs)); 37*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecSetUp(rhs)); 38c4762a1bSJed Brown 39c4762a1bSJed Brown rhsval = 0.0; 40c4762a1bSJed Brown for (i=0; i<m; i++) { 41c4762a1bSJed Brown for (j=nlocal*rank; j<nlocal*(rank+1); j++) { 42c4762a1bSJed Brown a = a0; 43c4762a1bSJed Brown for (b=0; b<bs; b++) { 44c4762a1bSJed Brown /* let's start with a 5-point stencil diffusion term */ 45c4762a1bSJed Brown v = -1.0; Ii = (j + n*i)*bs + b; 46*5f80ce2aSJacob Faibussowitsch if (i>0) {J = Ii - n*bs; CHKERRQ(MatSetValues(A,1,&Ii,1,&J,&v,ADD_VALUES));} 47*5f80ce2aSJacob Faibussowitsch if (i<m-1) {J = Ii + n*bs; CHKERRQ(MatSetValues(A,1,&Ii,1,&J,&v,ADD_VALUES));} 48*5f80ce2aSJacob Faibussowitsch if (j>0) {J = Ii - 1*bs; CHKERRQ(MatSetValues(A,1,&Ii,1,&J,&v,ADD_VALUES));} 49*5f80ce2aSJacob Faibussowitsch if (j<n-1) {J = Ii + 1*bs; CHKERRQ(MatSetValues(A,1,&Ii,1,&J,&v,ADD_VALUES));} 50*5f80ce2aSJacob Faibussowitsch v = 4.0; CHKERRQ(MatSetValues(A,1,&Ii,1,&Ii,&v,ADD_VALUES)); 51c4762a1bSJed Brown if (upwind) { 52c4762a1bSJed Brown /* now add a 2nd order upwind advection term to add a little asymmetry */ 53c4762a1bSJed Brown if (j>2) { 54c4762a1bSJed Brown J = Ii-2*bs; v2 = 0.5*a; v1 = -2.0*a; v0 = 1.5*a; 55*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValues(A,1,&Ii,1,&J,&v2,ADD_VALUES)); 56c4762a1bSJed Brown } else { 57c4762a1bSJed Brown /* fall back to 1st order upwind */ 58c4762a1bSJed Brown v1 = -1.0*a; v0 = 1.0*a; 59c4762a1bSJed Brown }; 60*5f80ce2aSJacob Faibussowitsch if (j>1) {J = Ii-1*bs; CHKERRQ(MatSetValues(A,1,&Ii,1,&J,&v1,ADD_VALUES));} 61*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValues(A,1,&Ii,1,&Ii,&v0,ADD_VALUES)); 62c4762a1bSJed Brown a /= 10.; /* use a different velocity for the next component */ 63c4762a1bSJed Brown /* add a coupling to the previous and next components */ 64c4762a1bSJed Brown v = 0.5; 65*5f80ce2aSJacob Faibussowitsch if (b>0) {J = Ii - 1; CHKERRQ(MatSetValues(A,1,&Ii,1,&J,&v,ADD_VALUES));} 66*5f80ce2aSJacob Faibussowitsch if (b<bs-1) {J = Ii + 1; CHKERRQ(MatSetValues(A,1,&Ii,1,&J,&v,ADD_VALUES));} 67c4762a1bSJed Brown } 68c4762a1bSJed Brown /* make up some rhs */ 69*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecSetValue(rhs, Ii, rhsval, INSERT_VALUES)); 70c4762a1bSJed Brown rhsval += 1.0; 71c4762a1bSJed Brown } 72c4762a1bSJed Brown } 73c4762a1bSJed Brown } 74*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY)); 75*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY)); 76c4762a1bSJed Brown 77c4762a1bSJed Brown if (convert) { /* Test different Mat implementations */ 78c4762a1bSJed Brown Mat B; 79c4762a1bSJed Brown 80*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatConvert(A,convname,MAT_INITIAL_MATRIX,&B)); 81*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&A)); 82c4762a1bSJed Brown A = B; 83c4762a1bSJed Brown } 84c4762a1bSJed Brown 85*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecAssemblyBegin(rhs)); 86*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecAssemblyEnd(rhs)); 87c4762a1bSJed Brown /* set rhs to zero to simplify */ 88c4762a1bSJed Brown if (zerorhs) { 89*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecZeroEntries(rhs)); 90c4762a1bSJed Brown } 91c4762a1bSJed Brown 92c4762a1bSJed Brown if (nonlocalBC) { 93c4762a1bSJed Brown /*version where boundary conditions are set by processes that don't necessarily own the nodes */ 94dd400576SPatrick Sanan if (rank == 0) { 95c4762a1bSJed Brown nboundary_nodes = size>m ? nlocal : m-size+nlocal; 96*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc1(nboundary_nodes,&boundary_nodes)); 97c4762a1bSJed Brown k = 0; 98c4762a1bSJed Brown for (i=size; i<m; i++,k++) {boundary_nodes[k] = n*i;}; 99c4762a1bSJed Brown } else if (rank < m) { 100c4762a1bSJed Brown nboundary_nodes = nlocal+1; 101*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc1(nboundary_nodes,&boundary_nodes)); 102c4762a1bSJed Brown boundary_nodes[0] = rank*n; 103c4762a1bSJed Brown k = 1; 104c4762a1bSJed Brown } else { 105c4762a1bSJed Brown nboundary_nodes = nlocal; 106*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc1(nboundary_nodes,&boundary_nodes)); 107c4762a1bSJed Brown k = 0; 108c4762a1bSJed Brown } 109c4762a1bSJed Brown for (j=nlocal*rank; j<nlocal*(rank+1); j++,k++) {boundary_nodes[k] = j;}; 110c4762a1bSJed Brown } else { 111c4762a1bSJed Brown /*version where boundary conditions are set by the node owners only */ 112*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc1(m*n,&boundary_nodes)); 113c4762a1bSJed Brown k=0; 114c4762a1bSJed Brown for (j=0; j<n; j++) { 115c4762a1bSJed Brown Ii = j; 116c4762a1bSJed Brown if (Ii>=rank*m*nlocal && Ii<(rank+1)*m*nlocal) boundary_nodes[k++] = Ii; 117c4762a1bSJed Brown } 118c4762a1bSJed Brown for (i=1; i<m; i++) { 119c4762a1bSJed Brown Ii = n*i; 120c4762a1bSJed Brown if (Ii>=rank*m*nlocal && Ii<(rank+1)*m*nlocal) boundary_nodes[k++] = Ii; 121c4762a1bSJed Brown } 122c4762a1bSJed Brown nboundary_nodes = k; 123c4762a1bSJed Brown } 124c4762a1bSJed Brown 125*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecDuplicate(rhs, &x)); 126*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecZeroEntries(x)); 127*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc2(nboundary_nodes*bs,&boundary_indices,nboundary_nodes*bs,&boundary_values)); 128c4762a1bSJed Brown for (k=0; k<nboundary_nodes; k++) { 129c4762a1bSJed Brown Ii = boundary_nodes[k]*bs; 130c4762a1bSJed Brown v = 1.0*boundary_nodes[k]; 131c4762a1bSJed Brown for (b=0; b<bs; b++, Ii++) { 132c4762a1bSJed Brown boundary_indices[k*bs+b] = Ii; 133c4762a1bSJed Brown boundary_values[k*bs+b] = v; 134*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscSynchronizedPrintf(PETSC_COMM_WORLD, "%d %" PetscInt_FMT " %f\n", rank, Ii, (double)PetscRealPart(v))); 135c4762a1bSJed Brown v += 0.1; 136c4762a1bSJed Brown } 137c4762a1bSJed Brown } 138*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscSynchronizedFlush(PETSC_COMM_WORLD, NULL)); 139*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecSetValues(x, nboundary_nodes*bs, boundary_indices, boundary_values, INSERT_VALUES)); 140*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecAssemblyBegin(x)); 141*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecAssemblyEnd(x)); 142c4762a1bSJed Brown 143c4762a1bSJed Brown /* We can check the rhs returned by MatZeroColumns by computing y=rhs-A*x and overwriting the boundary entries with boundary values */ 144*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecDuplicate(x, &y)); 145*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatMult(A, x, y)); 146*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecAYPX(y, -1.0, rhs)); 147c4762a1bSJed Brown for (k=0; k<nboundary_nodes*bs; k++) boundary_values[k] *= diag; 148*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecSetValues(y, nboundary_nodes*bs, boundary_indices, boundary_values, INSERT_VALUES)); 149*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecAssemblyBegin(y)); 150*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecAssemblyEnd(y)); 151c4762a1bSJed Brown 152*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_WORLD, "*** Matrix A and vector x:\n")); 153*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatView(A, PETSC_VIEWER_STDOUT_WORLD)); 154*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecView(x,PETSC_VIEWER_STDOUT_WORLD)); 155c4762a1bSJed Brown 156*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatZeroRowsColumns(A, nboundary_nodes*bs, boundary_indices, diag, x, rhs)); 157*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_WORLD, "*** Vector rhs returned by MatZeroRowsColumns\n")); 158*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecView(rhs,PETSC_VIEWER_STDOUT_WORLD)); 159*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecAXPY(y, -1.0, rhs)); 160*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecNorm(y, NORM_INFINITY, &norm)); 161c4762a1bSJed Brown if (norm > 1.0e-10) { 162*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_WORLD, "*** Difference between rhs and y, inf-norm: %f\n", (double)norm)); 163*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecView(y,PETSC_VIEWER_STDOUT_WORLD)); 164c4762a1bSJed Brown SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_PLIB, "Bug in MatZeroRowsColumns"); 165c4762a1bSJed Brown } 166c4762a1bSJed Brown 167*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree(boundary_nodes)); 168*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree2(boundary_indices,boundary_values)); 169*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&x)); 170*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&y)); 171*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&rhs)); 172*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&A)); 173c4762a1bSJed Brown 174c4762a1bSJed Brown ierr = PetscFinalize(); 175c4762a1bSJed Brown return ierr; 176c4762a1bSJed Brown } 177c4762a1bSJed Brown 178c4762a1bSJed Brown /*TEST 179c4762a1bSJed Brown 180c4762a1bSJed Brown test: 1819fd945daSStefano Zampini diff_args: -j 182c4762a1bSJed Brown suffix: 0 183c4762a1bSJed Brown 184c4762a1bSJed Brown test: 1859fd945daSStefano Zampini diff_args: -j 186c4762a1bSJed Brown suffix: 1 187c4762a1bSJed Brown nsize: 2 188c4762a1bSJed Brown 189c4762a1bSJed Brown test: 1909fd945daSStefano Zampini diff_args: -j 191c4762a1bSJed Brown suffix: 10 192c4762a1bSJed Brown nsize: 2 193c4762a1bSJed Brown args: -bs 2 -nonlocal_bc 194c4762a1bSJed Brown 195c4762a1bSJed Brown test: 1969fd945daSStefano Zampini diff_args: -j 197c4762a1bSJed Brown suffix: 11 198c4762a1bSJed Brown nsize: 7 199c4762a1bSJed Brown args: -bs 2 -nonlocal_bc 200c4762a1bSJed Brown 201c4762a1bSJed Brown test: 2029fd945daSStefano Zampini diff_args: -j 203c4762a1bSJed Brown suffix: 12 204c4762a1bSJed Brown args: -bs 2 -nonlocal_bc -mat_type baij 205c4762a1bSJed Brown 206c4762a1bSJed Brown test: 2079fd945daSStefano Zampini diff_args: -j 208c4762a1bSJed Brown suffix: 13 209c4762a1bSJed Brown nsize: 2 210c4762a1bSJed Brown args: -bs 2 -nonlocal_bc -mat_type baij 211c4762a1bSJed Brown 212c4762a1bSJed Brown test: 2139fd945daSStefano Zampini diff_args: -j 214c4762a1bSJed Brown suffix: 14 215c4762a1bSJed Brown nsize: 7 216c4762a1bSJed Brown args: -bs 2 -nonlocal_bc -mat_type baij 217c4762a1bSJed Brown 218c4762a1bSJed Brown test: 2199fd945daSStefano Zampini diff_args: -j 220c4762a1bSJed Brown suffix: 2 221c4762a1bSJed Brown nsize: 7 222c4762a1bSJed Brown 223c4762a1bSJed Brown test: 2249fd945daSStefano Zampini diff_args: -j 225c4762a1bSJed Brown suffix: 3 226c4762a1bSJed Brown args: -mat_type baij 227c4762a1bSJed Brown 228c4762a1bSJed Brown test: 2299fd945daSStefano Zampini diff_args: -j 230c4762a1bSJed Brown suffix: 4 231c4762a1bSJed Brown nsize: 2 232c4762a1bSJed Brown args: -mat_type baij 233c4762a1bSJed Brown 234c4762a1bSJed Brown test: 2359fd945daSStefano Zampini diff_args: -j 236c4762a1bSJed Brown suffix: 5 237c4762a1bSJed Brown nsize: 7 238c4762a1bSJed Brown args: -mat_type baij 239c4762a1bSJed Brown 240c4762a1bSJed Brown test: 2419fd945daSStefano Zampini diff_args: -j 242c4762a1bSJed Brown suffix: 6 243c4762a1bSJed Brown args: -bs 2 -mat_type baij 244c4762a1bSJed Brown 245c4762a1bSJed Brown test: 2469fd945daSStefano Zampini diff_args: -j 247c4762a1bSJed Brown suffix: 7 248c4762a1bSJed Brown nsize: 2 249c4762a1bSJed Brown args: -bs 2 -mat_type baij 250c4762a1bSJed Brown 251c4762a1bSJed Brown test: 2529fd945daSStefano Zampini diff_args: -j 253c4762a1bSJed Brown suffix: 8 254c4762a1bSJed Brown nsize: 7 255c4762a1bSJed Brown args: -bs 2 -mat_type baij 256c4762a1bSJed Brown 257c4762a1bSJed Brown test: 2589fd945daSStefano Zampini diff_args: -j 259c4762a1bSJed Brown suffix: 9 260c4762a1bSJed Brown args: -bs 2 -nonlocal_bc 261c4762a1bSJed Brown 262c4762a1bSJed Brown test: 2639fd945daSStefano Zampini diff_args: -j 264c4762a1bSJed Brown suffix: 15 265c4762a1bSJed Brown args: -bs 2 -nonlocal_bc -convname shell 266c4762a1bSJed Brown 267c4762a1bSJed Brown test: 2689fd945daSStefano Zampini diff_args: -j 269c4762a1bSJed Brown suffix: 16 270c4762a1bSJed Brown nsize: 2 271c4762a1bSJed Brown args: -bs 2 -nonlocal_bc -convname shell 272c4762a1bSJed Brown 273c4762a1bSJed Brown test: 2749fd945daSStefano Zampini diff_args: -j 275c4762a1bSJed Brown suffix: 17 276c4762a1bSJed Brown args: -bs 2 -nonlocal_bc -convname dense 277c4762a1bSJed Brown 278c4762a1bSJed Brown testset: 2799fd945daSStefano Zampini diff_args: -j 280c4762a1bSJed Brown suffix: full 281c4762a1bSJed Brown nsize: {{1 3}separate output} 282c4762a1bSJed Brown args: -diag {{0.12 -0.13}separate output} -convname {{aij shell baij}separate output} -zerorhs 0 28338a8e8c1SStefano Zampini 28438a8e8c1SStefano Zampini test: 2859fd945daSStefano Zampini diff_args: -j 28638a8e8c1SStefano Zampini requires: cuda 28738a8e8c1SStefano Zampini suffix: cusparse_1 28838a8e8c1SStefano Zampini nsize: 1 28938a8e8c1SStefano Zampini args: -diag {{0.12 -0.13}separate output} -convname {{seqaijcusparse mpiaijcusparse}separate output} -zerorhs 0 -mat_type {{seqaijcusparse mpiaijcusparse}separate output} 29038a8e8c1SStefano Zampini 29138a8e8c1SStefano Zampini test: 2929fd945daSStefano Zampini diff_args: -j 29338a8e8c1SStefano Zampini requires: cuda 29438a8e8c1SStefano Zampini suffix: cusparse_3 29538a8e8c1SStefano Zampini nsize: 3 29638a8e8c1SStefano Zampini args: -diag {{0.12 -0.13}separate output} -convname mpiaijcusparse -zerorhs 0 -mat_type mpiaijcusparse 29738a8e8c1SStefano Zampini 298c4762a1bSJed Brown TEST*/ 299