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 PetscScalar v,v0,v1,v2,a0=0.1,a,rhsval, *boundary_values,diag = 1.0; 14c4762a1bSJed Brown PetscReal norm; 15c4762a1bSJed Brown char convname[64]; 16c4762a1bSJed Brown PetscBool upwind = PETSC_FALSE, nonlocalBC = PETSC_FALSE, zerorhs = PETSC_TRUE, convert = PETSC_FALSE; 17c4762a1bSJed Brown 18*327415f7SBarry Smith PetscFunctionBeginUser; 199566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc,&args,(char*)0,help)); 209566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&rank)); 219566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size)); 22c4762a1bSJed Brown n = nlocal*size; 23c4762a1bSJed Brown 249566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL,NULL, "-bs", &bs, NULL)); 259566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL,NULL, "-nonlocal_bc", &nonlocalBC, NULL)); 269566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetScalar(NULL,NULL, "-diag", &diag, NULL)); 279566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetString(NULL,NULL,"-convname",convname,sizeof(convname),&convert)); 289566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL,NULL, "-zerorhs", &zerorhs, NULL)); 29c4762a1bSJed Brown 309566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD,&A)); 319566063dSJacob Faibussowitsch PetscCall(MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,m*n*bs,m*n*bs)); 329566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(A)); 339566063dSJacob Faibussowitsch PetscCall(MatSetUp(A)); 34c4762a1bSJed Brown 359566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(A, NULL, &rhs)); 369566063dSJacob Faibussowitsch PetscCall(VecSetFromOptions(rhs)); 379566063dSJacob Faibussowitsch PetscCall(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; 469566063dSJacob Faibussowitsch if (i>0) {J = Ii - n*bs; PetscCall(MatSetValues(A,1,&Ii,1,&J,&v,ADD_VALUES));} 479566063dSJacob Faibussowitsch if (i<m-1) {J = Ii + n*bs; PetscCall(MatSetValues(A,1,&Ii,1,&J,&v,ADD_VALUES));} 489566063dSJacob Faibussowitsch if (j>0) {J = Ii - 1*bs; PetscCall(MatSetValues(A,1,&Ii,1,&J,&v,ADD_VALUES));} 499566063dSJacob Faibussowitsch if (j<n-1) {J = Ii + 1*bs; PetscCall(MatSetValues(A,1,&Ii,1,&J,&v,ADD_VALUES));} 509566063dSJacob Faibussowitsch v = 4.0; PetscCall(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; 559566063dSJacob Faibussowitsch PetscCall(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 }; 609566063dSJacob Faibussowitsch if (j>1) {J = Ii-1*bs; PetscCall(MatSetValues(A,1,&Ii,1,&J,&v1,ADD_VALUES));} 619566063dSJacob Faibussowitsch PetscCall(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; 659566063dSJacob Faibussowitsch if (b>0) {J = Ii - 1; PetscCall(MatSetValues(A,1,&Ii,1,&J,&v,ADD_VALUES));} 669566063dSJacob Faibussowitsch if (b<bs-1) {J = Ii + 1; PetscCall(MatSetValues(A,1,&Ii,1,&J,&v,ADD_VALUES));} 67c4762a1bSJed Brown } 68c4762a1bSJed Brown /* make up some rhs */ 699566063dSJacob Faibussowitsch PetscCall(VecSetValue(rhs, Ii, rhsval, INSERT_VALUES)); 70c4762a1bSJed Brown rhsval += 1.0; 71c4762a1bSJed Brown } 72c4762a1bSJed Brown } 73c4762a1bSJed Brown } 749566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY)); 759566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY)); 76c4762a1bSJed Brown 77c4762a1bSJed Brown if (convert) { /* Test different Mat implementations */ 78c4762a1bSJed Brown Mat B; 79c4762a1bSJed Brown 809566063dSJacob Faibussowitsch PetscCall(MatConvert(A,convname,MAT_INITIAL_MATRIX,&B)); 819566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A)); 82c4762a1bSJed Brown A = B; 83c4762a1bSJed Brown } 84c4762a1bSJed Brown 859566063dSJacob Faibussowitsch PetscCall(VecAssemblyBegin(rhs)); 869566063dSJacob Faibussowitsch PetscCall(VecAssemblyEnd(rhs)); 87c4762a1bSJed Brown /* set rhs to zero to simplify */ 881baa6e33SBarry Smith if (zerorhs) PetscCall(VecZeroEntries(rhs)); 89c4762a1bSJed Brown 90c4762a1bSJed Brown if (nonlocalBC) { 91c4762a1bSJed Brown /*version where boundary conditions are set by processes that don't necessarily own the nodes */ 92dd400576SPatrick Sanan if (rank == 0) { 93c4762a1bSJed Brown nboundary_nodes = size>m ? nlocal : m-size+nlocal; 949566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nboundary_nodes,&boundary_nodes)); 95c4762a1bSJed Brown k = 0; 96c4762a1bSJed Brown for (i=size; i<m; i++,k++) {boundary_nodes[k] = n*i;}; 97c4762a1bSJed Brown } else if (rank < m) { 98c4762a1bSJed Brown nboundary_nodes = nlocal+1; 999566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nboundary_nodes,&boundary_nodes)); 100c4762a1bSJed Brown boundary_nodes[0] = rank*n; 101c4762a1bSJed Brown k = 1; 102c4762a1bSJed Brown } else { 103c4762a1bSJed Brown nboundary_nodes = nlocal; 1049566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nboundary_nodes,&boundary_nodes)); 105c4762a1bSJed Brown k = 0; 106c4762a1bSJed Brown } 107c4762a1bSJed Brown for (j=nlocal*rank; j<nlocal*(rank+1); j++,k++) {boundary_nodes[k] = j;}; 108c4762a1bSJed Brown } else { 109c4762a1bSJed Brown /*version where boundary conditions are set by the node owners only */ 1109566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(m*n,&boundary_nodes)); 111c4762a1bSJed Brown k=0; 112c4762a1bSJed Brown for (j=0; j<n; j++) { 113c4762a1bSJed Brown Ii = j; 114c4762a1bSJed Brown if (Ii>=rank*m*nlocal && Ii<(rank+1)*m*nlocal) boundary_nodes[k++] = Ii; 115c4762a1bSJed Brown } 116c4762a1bSJed Brown for (i=1; i<m; i++) { 117c4762a1bSJed Brown Ii = n*i; 118c4762a1bSJed Brown if (Ii>=rank*m*nlocal && Ii<(rank+1)*m*nlocal) boundary_nodes[k++] = Ii; 119c4762a1bSJed Brown } 120c4762a1bSJed Brown nboundary_nodes = k; 121c4762a1bSJed Brown } 122c4762a1bSJed Brown 1239566063dSJacob Faibussowitsch PetscCall(VecDuplicate(rhs, &x)); 1249566063dSJacob Faibussowitsch PetscCall(VecZeroEntries(x)); 1259566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(nboundary_nodes*bs,&boundary_indices,nboundary_nodes*bs,&boundary_values)); 126c4762a1bSJed Brown for (k=0; k<nboundary_nodes; k++) { 127c4762a1bSJed Brown Ii = boundary_nodes[k]*bs; 128c4762a1bSJed Brown v = 1.0*boundary_nodes[k]; 129c4762a1bSJed Brown for (b=0; b<bs; b++, Ii++) { 130c4762a1bSJed Brown boundary_indices[k*bs+b] = Ii; 131c4762a1bSJed Brown boundary_values[k*bs+b] = v; 1329566063dSJacob Faibussowitsch PetscCall(PetscSynchronizedPrintf(PETSC_COMM_WORLD, "%d %" PetscInt_FMT " %f\n", rank, Ii, (double)PetscRealPart(v))); 133c4762a1bSJed Brown v += 0.1; 134c4762a1bSJed Brown } 135c4762a1bSJed Brown } 1369566063dSJacob Faibussowitsch PetscCall(PetscSynchronizedFlush(PETSC_COMM_WORLD, NULL)); 1379566063dSJacob Faibussowitsch PetscCall(VecSetValues(x, nboundary_nodes*bs, boundary_indices, boundary_values, INSERT_VALUES)); 1389566063dSJacob Faibussowitsch PetscCall(VecAssemblyBegin(x)); 1399566063dSJacob Faibussowitsch PetscCall(VecAssemblyEnd(x)); 140c4762a1bSJed Brown 141c4762a1bSJed Brown /* We can check the rhs returned by MatZeroColumns by computing y=rhs-A*x and overwriting the boundary entries with boundary values */ 1429566063dSJacob Faibussowitsch PetscCall(VecDuplicate(x, &y)); 1439566063dSJacob Faibussowitsch PetscCall(MatMult(A, x, y)); 1449566063dSJacob Faibussowitsch PetscCall(VecAYPX(y, -1.0, rhs)); 145c4762a1bSJed Brown for (k=0; k<nboundary_nodes*bs; k++) boundary_values[k] *= diag; 1469566063dSJacob Faibussowitsch PetscCall(VecSetValues(y, nboundary_nodes*bs, boundary_indices, boundary_values, INSERT_VALUES)); 1479566063dSJacob Faibussowitsch PetscCall(VecAssemblyBegin(y)); 1489566063dSJacob Faibussowitsch PetscCall(VecAssemblyEnd(y)); 149c4762a1bSJed Brown 1509566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "*** Matrix A and vector x:\n")); 1519566063dSJacob Faibussowitsch PetscCall(MatView(A, PETSC_VIEWER_STDOUT_WORLD)); 1529566063dSJacob Faibussowitsch PetscCall(VecView(x,PETSC_VIEWER_STDOUT_WORLD)); 153c4762a1bSJed Brown 1549566063dSJacob Faibussowitsch PetscCall(MatZeroRowsColumns(A, nboundary_nodes*bs, boundary_indices, diag, x, rhs)); 1559566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "*** Vector rhs returned by MatZeroRowsColumns\n")); 1569566063dSJacob Faibussowitsch PetscCall(VecView(rhs,PETSC_VIEWER_STDOUT_WORLD)); 1579566063dSJacob Faibussowitsch PetscCall(VecAXPY(y, -1.0, rhs)); 1589566063dSJacob Faibussowitsch PetscCall(VecNorm(y, NORM_INFINITY, &norm)); 159c4762a1bSJed Brown if (norm > 1.0e-10) { 1609566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "*** Difference between rhs and y, inf-norm: %f\n", (double)norm)); 1619566063dSJacob Faibussowitsch PetscCall(VecView(y,PETSC_VIEWER_STDOUT_WORLD)); 162c4762a1bSJed Brown SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_PLIB, "Bug in MatZeroRowsColumns"); 163c4762a1bSJed Brown } 164c4762a1bSJed Brown 1659566063dSJacob Faibussowitsch PetscCall(PetscFree(boundary_nodes)); 1669566063dSJacob Faibussowitsch PetscCall(PetscFree2(boundary_indices,boundary_values)); 1679566063dSJacob Faibussowitsch PetscCall(VecDestroy(&x)); 1689566063dSJacob Faibussowitsch PetscCall(VecDestroy(&y)); 1699566063dSJacob Faibussowitsch PetscCall(VecDestroy(&rhs)); 1709566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A)); 171c4762a1bSJed Brown 1729566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 173b122ec5aSJacob Faibussowitsch return 0; 174c4762a1bSJed Brown } 175c4762a1bSJed Brown 176c4762a1bSJed Brown /*TEST 177c4762a1bSJed Brown 178c4762a1bSJed Brown test: 1799fd945daSStefano Zampini diff_args: -j 180c4762a1bSJed Brown suffix: 0 181c4762a1bSJed Brown 182c4762a1bSJed Brown test: 1839fd945daSStefano Zampini diff_args: -j 184c4762a1bSJed Brown suffix: 1 185c4762a1bSJed Brown nsize: 2 186c4762a1bSJed Brown 187c4762a1bSJed Brown test: 1889fd945daSStefano Zampini diff_args: -j 189c4762a1bSJed Brown suffix: 10 190c4762a1bSJed Brown nsize: 2 191c4762a1bSJed Brown args: -bs 2 -nonlocal_bc 192c4762a1bSJed Brown 193c4762a1bSJed Brown test: 1949fd945daSStefano Zampini diff_args: -j 195c4762a1bSJed Brown suffix: 11 196c4762a1bSJed Brown nsize: 7 197c4762a1bSJed Brown args: -bs 2 -nonlocal_bc 198c4762a1bSJed Brown 199c4762a1bSJed Brown test: 2009fd945daSStefano Zampini diff_args: -j 201c4762a1bSJed Brown suffix: 12 202c4762a1bSJed Brown args: -bs 2 -nonlocal_bc -mat_type baij 203c4762a1bSJed Brown 204c4762a1bSJed Brown test: 2059fd945daSStefano Zampini diff_args: -j 206c4762a1bSJed Brown suffix: 13 207c4762a1bSJed Brown nsize: 2 208c4762a1bSJed Brown args: -bs 2 -nonlocal_bc -mat_type baij 209c4762a1bSJed Brown 210c4762a1bSJed Brown test: 2119fd945daSStefano Zampini diff_args: -j 212c4762a1bSJed Brown suffix: 14 213c4762a1bSJed Brown nsize: 7 214c4762a1bSJed Brown args: -bs 2 -nonlocal_bc -mat_type baij 215c4762a1bSJed Brown 216c4762a1bSJed Brown test: 2179fd945daSStefano Zampini diff_args: -j 218c4762a1bSJed Brown suffix: 2 219c4762a1bSJed Brown nsize: 7 220c4762a1bSJed Brown 221c4762a1bSJed Brown test: 2229fd945daSStefano Zampini diff_args: -j 223c4762a1bSJed Brown suffix: 3 224c4762a1bSJed Brown args: -mat_type baij 225c4762a1bSJed Brown 226c4762a1bSJed Brown test: 2279fd945daSStefano Zampini diff_args: -j 228c4762a1bSJed Brown suffix: 4 229c4762a1bSJed Brown nsize: 2 230c4762a1bSJed Brown args: -mat_type baij 231c4762a1bSJed Brown 232c4762a1bSJed Brown test: 2339fd945daSStefano Zampini diff_args: -j 234c4762a1bSJed Brown suffix: 5 235c4762a1bSJed Brown nsize: 7 236c4762a1bSJed Brown args: -mat_type baij 237c4762a1bSJed Brown 238c4762a1bSJed Brown test: 2399fd945daSStefano Zampini diff_args: -j 240c4762a1bSJed Brown suffix: 6 241c4762a1bSJed Brown args: -bs 2 -mat_type baij 242c4762a1bSJed Brown 243c4762a1bSJed Brown test: 2449fd945daSStefano Zampini diff_args: -j 245c4762a1bSJed Brown suffix: 7 246c4762a1bSJed Brown nsize: 2 247c4762a1bSJed Brown args: -bs 2 -mat_type baij 248c4762a1bSJed Brown 249c4762a1bSJed Brown test: 2509fd945daSStefano Zampini diff_args: -j 251c4762a1bSJed Brown suffix: 8 252c4762a1bSJed Brown nsize: 7 253c4762a1bSJed Brown args: -bs 2 -mat_type baij 254c4762a1bSJed Brown 255c4762a1bSJed Brown test: 2569fd945daSStefano Zampini diff_args: -j 257c4762a1bSJed Brown suffix: 9 258c4762a1bSJed Brown args: -bs 2 -nonlocal_bc 259c4762a1bSJed Brown 260c4762a1bSJed Brown test: 2619fd945daSStefano Zampini diff_args: -j 262c4762a1bSJed Brown suffix: 15 263c4762a1bSJed Brown args: -bs 2 -nonlocal_bc -convname shell 264c4762a1bSJed Brown 265c4762a1bSJed Brown test: 2669fd945daSStefano Zampini diff_args: -j 267c4762a1bSJed Brown suffix: 16 268c4762a1bSJed Brown nsize: 2 269c4762a1bSJed Brown args: -bs 2 -nonlocal_bc -convname shell 270c4762a1bSJed Brown 271c4762a1bSJed Brown test: 2729fd945daSStefano Zampini diff_args: -j 273c4762a1bSJed Brown suffix: 17 274c4762a1bSJed Brown args: -bs 2 -nonlocal_bc -convname dense 275c4762a1bSJed Brown 276c4762a1bSJed Brown testset: 2779fd945daSStefano Zampini diff_args: -j 278c4762a1bSJed Brown suffix: full 279c4762a1bSJed Brown nsize: {{1 3}separate output} 280c4762a1bSJed Brown args: -diag {{0.12 -0.13}separate output} -convname {{aij shell baij}separate output} -zerorhs 0 28138a8e8c1SStefano Zampini 28238a8e8c1SStefano Zampini test: 2839fd945daSStefano Zampini diff_args: -j 28438a8e8c1SStefano Zampini requires: cuda 28538a8e8c1SStefano Zampini suffix: cusparse_1 28638a8e8c1SStefano Zampini nsize: 1 28738a8e8c1SStefano Zampini args: -diag {{0.12 -0.13}separate output} -convname {{seqaijcusparse mpiaijcusparse}separate output} -zerorhs 0 -mat_type {{seqaijcusparse mpiaijcusparse}separate output} 28838a8e8c1SStefano Zampini 28938a8e8c1SStefano Zampini test: 2909fd945daSStefano Zampini diff_args: -j 29138a8e8c1SStefano Zampini requires: cuda 29238a8e8c1SStefano Zampini suffix: cusparse_3 29338a8e8c1SStefano Zampini nsize: 3 29438a8e8c1SStefano Zampini args: -diag {{0.12 -0.13}separate output} -convname mpiaijcusparse -zerorhs 0 -mat_type mpiaijcusparse 29538a8e8c1SStefano Zampini 296c4762a1bSJed Brown TEST*/ 297