1c4762a1bSJed Brown static char help[] = "Tests DMSLICED operations\n\n"; 2c4762a1bSJed Brown 3c4762a1bSJed Brown #include <petscdmsliced.h> 4c4762a1bSJed Brown 5c4762a1bSJed Brown int main(int argc,char *argv[]) 6c4762a1bSJed Brown { 7c4762a1bSJed Brown char mat_type[256] = MATAIJ; /* default matrix type */ 8c4762a1bSJed Brown PetscErrorCode ierr; 9c4762a1bSJed Brown MPI_Comm comm; 10c4762a1bSJed Brown PetscMPIInt rank,size; 11c4762a1bSJed Brown DM slice; 12c4762a1bSJed Brown PetscInt i,bs=1,N=5,n,m,rstart,ghosts[2],*d_nnz,*o_nnz,dfill[4]={1,0,0,1},ofill[4]={1,1,1,1}; 13c4762a1bSJed Brown PetscReal alpha =1,K=1,rho0=1,u0=0,sigma=0.2; 14c4762a1bSJed Brown PetscBool useblock=PETSC_TRUE; 15c4762a1bSJed Brown PetscScalar *xx; 16c4762a1bSJed Brown Mat A; 17c4762a1bSJed Brown Vec x,b,lf; 18c4762a1bSJed Brown 199566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc,&argv,0,help)); 20c4762a1bSJed Brown comm = PETSC_COMM_WORLD; 219566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(comm,&size)); 229566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(comm,&rank)); 23c4762a1bSJed Brown 249566063dSJacob Faibussowitsch ierr = PetscOptionsBegin(comm,0,"Options for DMSliced test",0);PetscCall(ierr); 25c4762a1bSJed Brown { 269566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-n","Global number of nodes","",N,&N,NULL)); 279566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-bs","Block size (1 or 2)","",bs,&bs,NULL)); 28c4762a1bSJed Brown if (bs != 1) { 29*08401ef6SPierre Jolivet PetscCheck(bs == 2,PETSC_COMM_WORLD,PETSC_ERR_SUP,"Block size must be 1 or 2"); 309566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-alpha","Inverse time step for wave operator","",alpha,&alpha,NULL)); 319566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-K","Bulk modulus of compressibility","",K,&K,NULL)); 329566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-rho0","Reference density","",rho0,&rho0,NULL)); 339566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-u0","Reference velocity","",u0,&u0,NULL)); 349566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-sigma","Width of Gaussian density perturbation","",sigma,&sigma,NULL)); 359566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-block","Use block matrix assembly","",useblock,&useblock,NULL)); 36c4762a1bSJed Brown } 379566063dSJacob Faibussowitsch PetscCall(PetscOptionsString("-sliced_mat_type","Matrix type to use (aij or baij)","",mat_type,mat_type,sizeof(mat_type),NULL)); 38c4762a1bSJed Brown } 399566063dSJacob Faibussowitsch ierr = PetscOptionsEnd();PetscCall(ierr); 40c4762a1bSJed Brown 41c4762a1bSJed Brown /* Split ownership, set up periodic grid in 1D */ 42c4762a1bSJed Brown n = PETSC_DECIDE; 439566063dSJacob Faibussowitsch PetscCall(PetscSplitOwnership(comm,&n,&N)); 44c4762a1bSJed Brown rstart = 0; 459566063dSJacob Faibussowitsch PetscCallMPI(MPI_Scan(&n,&rstart,1,MPIU_INT,MPI_SUM,comm)); 46c4762a1bSJed Brown rstart -= n; 47c4762a1bSJed Brown ghosts[0] = (N+rstart-1)%N; 48c4762a1bSJed Brown ghosts[1] = (rstart+n)%N; 49c4762a1bSJed Brown 509566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(n,&d_nnz,n,&o_nnz)); 51c4762a1bSJed Brown for (i=0; i<n; i++) { 52c4762a1bSJed Brown if (size > 1 && (i==0 || i==n-1)) { 53c4762a1bSJed Brown d_nnz[i] = 2; 54c4762a1bSJed Brown o_nnz[i] = 1; 55c4762a1bSJed Brown } else { 56c4762a1bSJed Brown d_nnz[i] = 3; 57c4762a1bSJed Brown o_nnz[i] = 0; 58c4762a1bSJed Brown } 59c4762a1bSJed Brown } 609566063dSJacob Faibussowitsch PetscCall(DMSlicedCreate(comm,bs,n,2,ghosts,d_nnz,o_nnz,&slice)); /* Currently does not copy X_nnz so we can't free them until after DMSlicedGetMatrix */ 61c4762a1bSJed Brown 629566063dSJacob Faibussowitsch if (!useblock) PetscCall(DMSlicedSetBlockFills(slice,dfill,ofill)); /* Irrelevant for baij formats */ 639566063dSJacob Faibussowitsch PetscCall(DMSetMatType(slice,mat_type)); 649566063dSJacob Faibussowitsch PetscCall(DMCreateMatrix(slice,&A)); 659566063dSJacob Faibussowitsch PetscCall(PetscFree2(d_nnz,o_nnz)); 669566063dSJacob Faibussowitsch PetscCall(MatSetOption(A,MAT_NEW_NONZERO_ALLOCATION_ERR,PETSC_TRUE)); 67c4762a1bSJed Brown 689566063dSJacob Faibussowitsch PetscCall(DMCreateGlobalVector(slice,&x)); 699566063dSJacob Faibussowitsch PetscCall(VecDuplicate(x,&b)); 70c4762a1bSJed Brown 719566063dSJacob Faibussowitsch PetscCall(VecGhostGetLocalForm(x,&lf)); 729566063dSJacob Faibussowitsch PetscCall(VecGetSize(lf,&m)); 73*08401ef6SPierre Jolivet PetscCheck(m == (n+2)*bs,PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"size of local form %D, expected %D",m,(n+2)*bs); 749566063dSJacob Faibussowitsch PetscCall(VecGetArray(lf,&xx)); 75c4762a1bSJed Brown for (i=0; i<n; i++) { 76c4762a1bSJed Brown PetscInt row[2],col[9],im,ip; 77c4762a1bSJed Brown PetscScalar v[12]; 78c4762a1bSJed Brown const PetscReal xref = 2.0*(rstart+i)/N - 1; /* [-1,1] */ 79c4762a1bSJed Brown const PetscReal h = 1.0/N; /* grid spacing */ 80c4762a1bSJed Brown im = (i==0) ? n : i-1; 81c4762a1bSJed Brown ip = (i==n-1) ? n+1 : i+1; 82c4762a1bSJed Brown switch (bs) { 83c4762a1bSJed Brown case 1: /* Laplacian with periodic boundaries */ 84c4762a1bSJed Brown col[0] = im; col[1] = i; col[2] = ip; 85c4762a1bSJed Brown v[0] = -h; v[1] = 2*h; v[2] = -h; 869566063dSJacob Faibussowitsch PetscCall(MatSetValuesLocal(A,1,&i,3,col,v,INSERT_VALUES)); 87c4762a1bSJed Brown xx[i] = PetscSinReal(xref*PETSC_PI); 88c4762a1bSJed Brown break; 89c4762a1bSJed Brown case 2: /* Linear acoustic wave operator in variables [rho, u], central differences, periodic, timestep 1/alpha */ 90c4762a1bSJed Brown v[0] = -0.5*u0; v[1] = -0.5*K; v[2] = alpha; v[3] = 0; v[4] = 0.5*u0; v[5] = 0.5*K; 91c4762a1bSJed Brown v[6] = -0.5/rho0; v[7] = -0.5*u0; v[8] = 0; v[9] = alpha; v[10] = 0.5/rho0; v[11] = 0.5*u0; 92c4762a1bSJed Brown if (useblock) { 93c4762a1bSJed Brown row[0] = i; col[0] = im; col[1] = i; col[2] = ip; 949566063dSJacob Faibussowitsch PetscCall(MatSetValuesBlockedLocal(A,1,row,3,col,v,INSERT_VALUES)); 95c4762a1bSJed Brown } else { 96c4762a1bSJed Brown row[0] = 2*i; row[1] = 2*i+1; 97c4762a1bSJed Brown col[0] = 2*im; col[1] = 2*im+1; col[2] = 2*i; col[3] = 2*ip; col[4] = 2*ip+1; 98c4762a1bSJed Brown v[3] = v[4]; v[4] = v[5]; /* pack values in first row */ 999566063dSJacob Faibussowitsch PetscCall(MatSetValuesLocal(A,1,row,5,col,v,INSERT_VALUES)); 100c4762a1bSJed Brown col[2] = 2*i+1; 101c4762a1bSJed Brown v[8] = v[9]; v[9] = v[10]; v[10] = v[11]; /* pack values in second row */ 1029566063dSJacob Faibussowitsch PetscCall(MatSetValuesLocal(A,1,row+1,5,col,v+6,INSERT_VALUES)); 103c4762a1bSJed Brown } 104c4762a1bSJed Brown /* Set current state (gaussian density perturbation) */ 105c4762a1bSJed Brown xx[2*i] = 0.2*PetscExpReal(-PetscSqr(xref)/(2*PetscSqr(sigma))); 106c4762a1bSJed Brown xx[2*i+1] = 0; 107c4762a1bSJed Brown break; 10898921bdaSJacob Faibussowitsch default: SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"not implemented for block size %D",bs); 109c4762a1bSJed Brown } 110c4762a1bSJed Brown } 1119566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(lf,&xx)); 1129566063dSJacob Faibussowitsch PetscCall(VecGhostRestoreLocalForm(x,&lf)); 1139566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY)); 1149566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY)); 115c4762a1bSJed Brown 1169566063dSJacob Faibussowitsch PetscCall(MatMult(A,x,b)); 1179566063dSJacob Faibussowitsch PetscCall(MatView(A,PETSC_VIEWER_STDOUT_WORLD)); 1189566063dSJacob Faibussowitsch PetscCall(VecView(x,PETSC_VIEWER_STDOUT_WORLD)); 1199566063dSJacob Faibussowitsch PetscCall(VecView(b,PETSC_VIEWER_STDOUT_WORLD)); 120c4762a1bSJed Brown 121c4762a1bSJed Brown /* Update the ghosted values, view the result on rank 0. */ 1229566063dSJacob Faibussowitsch PetscCall(VecGhostUpdateBegin(b,INSERT_VALUES,SCATTER_FORWARD)); 1239566063dSJacob Faibussowitsch PetscCall(VecGhostUpdateEnd(b,INSERT_VALUES,SCATTER_FORWARD)); 124dd400576SPatrick Sanan if (rank == 0) { 1259566063dSJacob Faibussowitsch PetscCall(VecGhostGetLocalForm(b,&lf)); 1269566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_SELF,"Local form of b on rank 0, last two nodes are ghost nodes\n")); 1279566063dSJacob Faibussowitsch PetscCall(VecView(lf,PETSC_VIEWER_STDOUT_SELF)); 1289566063dSJacob Faibussowitsch PetscCall(VecGhostRestoreLocalForm(b,&lf)); 129c4762a1bSJed Brown } 130c4762a1bSJed Brown 1319566063dSJacob Faibussowitsch PetscCall(DMDestroy(&slice)); 1329566063dSJacob Faibussowitsch PetscCall(VecDestroy(&x)); 1339566063dSJacob Faibussowitsch PetscCall(VecDestroy(&b)); 1349566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A)); 1359566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 136b122ec5aSJacob Faibussowitsch return 0; 137c4762a1bSJed Brown } 138c4762a1bSJed Brown 139c4762a1bSJed Brown /*TEST 140c4762a1bSJed Brown 141c4762a1bSJed Brown test: 142c4762a1bSJed Brown nsize: 2 143c4762a1bSJed Brown args: -bs 2 -block 0 -sliced_mat_type baij -alpha 10 -u0 0.1 144c4762a1bSJed Brown 145c4762a1bSJed Brown test: 146c4762a1bSJed Brown suffix: 2 147c4762a1bSJed Brown nsize: 2 148c4762a1bSJed Brown args: -bs 2 -block 1 -sliced_mat_type aij -alpha 10 -u0 0.1 149c4762a1bSJed Brown 150c4762a1bSJed Brown test: 151c4762a1bSJed Brown suffix: 3 152c4762a1bSJed Brown nsize: 2 153c4762a1bSJed Brown args: -bs 2 -block 0 -sliced_mat_type aij -alpha 10 -u0 0.1 154c4762a1bSJed Brown 155c4762a1bSJed Brown TEST*/ 156