1c4762a1bSJed Brown 2c4762a1bSJed Brown /* 3c4762a1bSJed Brown Laplacian in 3D. Use for testing MatSolve routines. 4c4762a1bSJed Brown Modeled by the partial differential equation 5c4762a1bSJed Brown 6c4762a1bSJed Brown - Laplacian u = 1,0 < x,y,z < 1, 7c4762a1bSJed Brown 8c4762a1bSJed Brown with boundary conditions 9c4762a1bSJed Brown u = 1 for x = 0, x = 1, y = 0, y = 1, z = 0, z = 1. 10c4762a1bSJed Brown */ 11c4762a1bSJed Brown 12c4762a1bSJed Brown static char help[] = "This example is for testing different MatSolve routines :MatSolve(), MatSolveAdd(), MatSolveTranspose(), MatSolveTransposeAdd(), and MatMatSolve().\n\ 13c4762a1bSJed Brown Example usage: ./ex129 -mat_type aij -dof 2\n\n"; 14c4762a1bSJed Brown 15c4762a1bSJed Brown #include <petscdm.h> 16c4762a1bSJed Brown #include <petscdmda.h> 17c4762a1bSJed Brown 18c4762a1bSJed Brown extern PetscErrorCode ComputeMatrix(DM,Mat); 19c4762a1bSJed Brown extern PetscErrorCode ComputeRHS(DM,Vec); 20c4762a1bSJed Brown extern PetscErrorCode ComputeRHSMatrix(PetscInt,PetscInt,Mat*); 21c4762a1bSJed Brown 22c4762a1bSJed Brown int main(int argc,char **args) 23c4762a1bSJed Brown { 24c4762a1bSJed Brown PetscMPIInt size; 25c4762a1bSJed Brown Vec x,b,y,b1; 26c4762a1bSJed Brown DM da; 27c4762a1bSJed Brown Mat A,F,RHS,X,C1; 28c4762a1bSJed Brown MatFactorInfo info; 29c4762a1bSJed Brown IS perm,iperm; 30c4762a1bSJed Brown PetscInt dof =1,M=8,m,n,nrhs; 31c4762a1bSJed Brown PetscScalar one = 1.0; 32c4762a1bSJed Brown PetscReal norm,tol = 1000*PETSC_MACHINE_EPSILON; 33c4762a1bSJed Brown PetscBool InplaceLU=PETSC_FALSE; 34c4762a1bSJed Brown 35*9566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc,&args,(char*)0,help)); 36*9566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size)); 372c71b3e2SJacob Faibussowitsch PetscCheckFalse(size != 1,PETSC_COMM_WORLD,PETSC_ERR_WRONG_MPI_SIZE,"This is a uniprocessor example only"); 38*9566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL,NULL,"-dof",&dof,NULL)); 39*9566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL,NULL,"-M",&M,NULL)); 40c4762a1bSJed Brown 41*9566063dSJacob Faibussowitsch PetscCall(DMDACreate(PETSC_COMM_WORLD,&da)); 42*9566063dSJacob Faibussowitsch PetscCall(DMSetDimension(da,3)); 43*9566063dSJacob Faibussowitsch PetscCall(DMDASetBoundaryType(da,DM_BOUNDARY_NONE,DM_BOUNDARY_NONE,DM_BOUNDARY_NONE)); 44*9566063dSJacob Faibussowitsch PetscCall(DMDASetStencilType(da,DMDA_STENCIL_STAR)); 45*9566063dSJacob Faibussowitsch PetscCall(DMDASetSizes(da,M,M,M)); 46*9566063dSJacob Faibussowitsch PetscCall(DMDASetNumProcs(da,PETSC_DECIDE,PETSC_DECIDE,PETSC_DECIDE)); 47*9566063dSJacob Faibussowitsch PetscCall(DMDASetDof(da,dof)); 48*9566063dSJacob Faibussowitsch PetscCall(DMDASetStencilWidth(da,1)); 49*9566063dSJacob Faibussowitsch PetscCall(DMDASetOwnershipRanges(da,NULL,NULL,NULL)); 50*9566063dSJacob Faibussowitsch PetscCall(DMSetMatType(da,MATBAIJ)); 51*9566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(da)); 52*9566063dSJacob Faibussowitsch PetscCall(DMSetUp(da)); 53c4762a1bSJed Brown 54*9566063dSJacob Faibussowitsch PetscCall(DMCreateGlobalVector(da,&x)); 55*9566063dSJacob Faibussowitsch PetscCall(DMCreateGlobalVector(da,&b)); 56*9566063dSJacob Faibussowitsch PetscCall(VecDuplicate(b,&y)); 57*9566063dSJacob Faibussowitsch PetscCall(ComputeRHS(da,b)); 58*9566063dSJacob Faibussowitsch PetscCall(VecSet(y,one)); 59*9566063dSJacob Faibussowitsch PetscCall(DMCreateMatrix(da,&A)); 60*9566063dSJacob Faibussowitsch PetscCall(ComputeMatrix(da,A)); 61*9566063dSJacob Faibussowitsch PetscCall(MatGetSize(A,&m,&n)); 62c4762a1bSJed Brown nrhs = 2; 63*9566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL,NULL,"-nrhs",&nrhs,NULL)); 64*9566063dSJacob Faibussowitsch PetscCall(ComputeRHSMatrix(m,nrhs,&RHS)); 65*9566063dSJacob Faibussowitsch PetscCall(MatDuplicate(RHS,MAT_DO_NOT_COPY_VALUES,&X)); 66c4762a1bSJed Brown 67*9566063dSJacob Faibussowitsch PetscCall(MatGetOrdering(A,MATORDERINGND,&perm,&iperm)); 68c4762a1bSJed Brown 69*9566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL,NULL,"-inplacelu",&InplaceLU,NULL)); 70*9566063dSJacob Faibussowitsch PetscCall(MatFactorInfoInitialize(&info)); 71c4762a1bSJed Brown if (!InplaceLU) { 72*9566063dSJacob Faibussowitsch PetscCall(MatGetFactor(A,MATSOLVERPETSC,MAT_FACTOR_LU,&F)); 73c4762a1bSJed Brown info.fill = 5.0; 74*9566063dSJacob Faibussowitsch PetscCall(MatLUFactorSymbolic(F,A,perm,iperm,&info)); 75*9566063dSJacob Faibussowitsch PetscCall(MatLUFactorNumeric(F,A,&info)); 76c4762a1bSJed Brown } else { /* Test inplace factorization */ 77*9566063dSJacob Faibussowitsch PetscCall(MatDuplicate(A,MAT_COPY_VALUES,&F)); 78*9566063dSJacob Faibussowitsch PetscCall(MatLUFactor(F,perm,iperm,&info)); 79c4762a1bSJed Brown } 80c4762a1bSJed Brown 81*9566063dSJacob Faibussowitsch PetscCall(VecDuplicate(y,&b1)); 82c4762a1bSJed Brown 83c4762a1bSJed Brown /* MatSolve */ 84*9566063dSJacob Faibussowitsch PetscCall(MatSolve(F,b,x)); 85*9566063dSJacob Faibussowitsch PetscCall(MatMult(A,x,b1)); 86*9566063dSJacob Faibussowitsch PetscCall(VecAXPY(b1,-1.0,b)); 87*9566063dSJacob Faibussowitsch PetscCall(VecNorm(b1,NORM_2,&norm)); 88c4762a1bSJed Brown if (norm > tol) { 89*9566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD,"MatSolve : Error of norm %g\n",(double)norm)); 90c4762a1bSJed Brown } 91c4762a1bSJed Brown 92c4762a1bSJed Brown /* MatSolveTranspose */ 93*9566063dSJacob Faibussowitsch PetscCall(MatSolveTranspose(F,b,x)); 94*9566063dSJacob Faibussowitsch PetscCall(MatMultTranspose(A,x,b1)); 95*9566063dSJacob Faibussowitsch PetscCall(VecAXPY(b1,-1.0,b)); 96*9566063dSJacob Faibussowitsch PetscCall(VecNorm(b1,NORM_2,&norm)); 97c4762a1bSJed Brown if (norm > tol) { 98*9566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD,"MatSolveTranspose : Error of norm %g\n",(double)norm)); 99c4762a1bSJed Brown } 100c4762a1bSJed Brown 101c4762a1bSJed Brown /* MatSolveAdd */ 102*9566063dSJacob Faibussowitsch PetscCall(MatSolveAdd(F,b,y,x)); 103*9566063dSJacob Faibussowitsch PetscCall(MatMult(A,y,b1)); 104*9566063dSJacob Faibussowitsch PetscCall(VecScale(b1,-1.0)); 105*9566063dSJacob Faibussowitsch PetscCall(MatMultAdd(A,x,b1,b1)); 106*9566063dSJacob Faibussowitsch PetscCall(VecAXPY(b1,-1.0,b)); 107*9566063dSJacob Faibussowitsch PetscCall(VecNorm(b1,NORM_2,&norm)); 108c4762a1bSJed Brown if (norm > tol) { 109*9566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD,"MatSolveAdd : Error of norm %g\n",(double)norm)); 110c4762a1bSJed Brown } 111c4762a1bSJed Brown 112c4762a1bSJed Brown /* MatSolveTransposeAdd */ 113*9566063dSJacob Faibussowitsch PetscCall(MatSolveTransposeAdd(F,b,y,x)); 114*9566063dSJacob Faibussowitsch PetscCall(MatMultTranspose(A,y,b1)); 115*9566063dSJacob Faibussowitsch PetscCall(VecScale(b1,-1.0)); 116*9566063dSJacob Faibussowitsch PetscCall(MatMultTransposeAdd(A,x,b1,b1)); 117*9566063dSJacob Faibussowitsch PetscCall(VecAXPY(b1,-1.0,b)); 118*9566063dSJacob Faibussowitsch PetscCall(VecNorm(b1,NORM_2,&norm)); 119c4762a1bSJed Brown if (norm > tol) { 120*9566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD,"MatSolveTransposeAdd : Error of norm %g\n",(double)norm)); 121c4762a1bSJed Brown } 122c4762a1bSJed Brown 123c4762a1bSJed Brown /* MatMatSolve */ 124*9566063dSJacob Faibussowitsch PetscCall(MatMatSolve(F,RHS,X)); 125*9566063dSJacob Faibussowitsch PetscCall(MatMatMult(A,X,MAT_INITIAL_MATRIX,2.0,&C1)); 126*9566063dSJacob Faibussowitsch PetscCall(MatAXPY(C1,-1.0,RHS,SAME_NONZERO_PATTERN)); 127*9566063dSJacob Faibussowitsch PetscCall(MatNorm(C1,NORM_FROBENIUS,&norm)); 128c4762a1bSJed Brown if (norm > tol) { 129*9566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD,"MatMatSolve : Error of norm %g\n",(double)norm)); 130c4762a1bSJed Brown } 131c4762a1bSJed Brown 132*9566063dSJacob Faibussowitsch PetscCall(VecDestroy(&x)); 133*9566063dSJacob Faibussowitsch PetscCall(VecDestroy(&b)); 134*9566063dSJacob Faibussowitsch PetscCall(VecDestroy(&b1)); 135*9566063dSJacob Faibussowitsch PetscCall(VecDestroy(&y)); 136*9566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A)); 137*9566063dSJacob Faibussowitsch PetscCall(MatDestroy(&F)); 138*9566063dSJacob Faibussowitsch PetscCall(MatDestroy(&RHS)); 139*9566063dSJacob Faibussowitsch PetscCall(MatDestroy(&C1)); 140*9566063dSJacob Faibussowitsch PetscCall(MatDestroy(&X)); 141*9566063dSJacob Faibussowitsch PetscCall(ISDestroy(&perm)); 142*9566063dSJacob Faibussowitsch PetscCall(ISDestroy(&iperm)); 143*9566063dSJacob Faibussowitsch PetscCall(DMDestroy(&da)); 144*9566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 145b122ec5aSJacob Faibussowitsch return 0; 146c4762a1bSJed Brown } 147c4762a1bSJed Brown 148c4762a1bSJed Brown PetscErrorCode ComputeRHS(DM da,Vec b) 149c4762a1bSJed Brown { 150c4762a1bSJed Brown PetscInt mx,my,mz; 151c4762a1bSJed Brown PetscScalar h; 152c4762a1bSJed Brown 153c4762a1bSJed Brown PetscFunctionBegin; 154*9566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(da,0,&mx,&my,&mz,0,0,0,0,0,0,0,0,0)); 155c4762a1bSJed Brown h = 1.0/((mx-1)*(my-1)*(mz-1)); 156*9566063dSJacob Faibussowitsch PetscCall(VecSet(b,h)); 157c4762a1bSJed Brown PetscFunctionReturn(0); 158c4762a1bSJed Brown } 159c4762a1bSJed Brown 160c4762a1bSJed Brown PetscErrorCode ComputeRHSMatrix(PetscInt m,PetscInt nrhs,Mat *C) 161c4762a1bSJed Brown { 162c4762a1bSJed Brown PetscRandom rand; 163c4762a1bSJed Brown Mat RHS; 164c4762a1bSJed Brown PetscScalar *array,rval; 165c4762a1bSJed Brown PetscInt i,k; 166c4762a1bSJed Brown 167c4762a1bSJed Brown PetscFunctionBegin; 168*9566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD,&RHS)); 169*9566063dSJacob Faibussowitsch PetscCall(MatSetSizes(RHS,m,PETSC_DECIDE,PETSC_DECIDE,nrhs)); 170*9566063dSJacob Faibussowitsch PetscCall(MatSetType(RHS,MATSEQDENSE)); 171*9566063dSJacob Faibussowitsch PetscCall(MatSetUp(RHS)); 172c4762a1bSJed Brown 173*9566063dSJacob Faibussowitsch PetscCall(PetscRandomCreate(PETSC_COMM_WORLD,&rand)); 174*9566063dSJacob Faibussowitsch PetscCall(PetscRandomSetFromOptions(rand)); 175*9566063dSJacob Faibussowitsch PetscCall(MatDenseGetArray(RHS,&array)); 176c4762a1bSJed Brown for (i=0; i<m; i++) { 177*9566063dSJacob Faibussowitsch PetscCall(PetscRandomGetValue(rand,&rval)); 178c4762a1bSJed Brown array[i] = rval; 179c4762a1bSJed Brown } 180c4762a1bSJed Brown if (nrhs > 1) { 181c4762a1bSJed Brown for (k=1; k<nrhs; k++) { 182c4762a1bSJed Brown for (i=0; i<m; i++) { 183c4762a1bSJed Brown array[m*k+i] = array[i]; 184c4762a1bSJed Brown } 185c4762a1bSJed Brown } 186c4762a1bSJed Brown } 187*9566063dSJacob Faibussowitsch PetscCall(MatDenseRestoreArray(RHS,&array)); 188*9566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(RHS,MAT_FINAL_ASSEMBLY)); 189*9566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(RHS,MAT_FINAL_ASSEMBLY)); 190c4762a1bSJed Brown *C = RHS; 191*9566063dSJacob Faibussowitsch PetscCall(PetscRandomDestroy(&rand)); 192c4762a1bSJed Brown PetscFunctionReturn(0); 193c4762a1bSJed Brown } 194c4762a1bSJed Brown 195c4762a1bSJed Brown PetscErrorCode ComputeMatrix(DM da,Mat B) 196c4762a1bSJed Brown { 197c4762a1bSJed Brown PetscInt i,j,k,mx,my,mz,xm,ym,zm,xs,ys,zs,dof,k1,k2,k3; 198c4762a1bSJed Brown PetscScalar *v,*v_neighbor,Hx,Hy,Hz,HxHydHz,HyHzdHx,HxHzdHy,r1,r2; 199c4762a1bSJed Brown MatStencil row,col; 200c4762a1bSJed Brown PetscRandom rand; 201c4762a1bSJed Brown 202c4762a1bSJed Brown PetscFunctionBegin; 203*9566063dSJacob Faibussowitsch PetscCall(PetscRandomCreate(PETSC_COMM_WORLD,&rand)); 204*9566063dSJacob Faibussowitsch PetscCall(PetscRandomSetSeed(rand,1)); 205*9566063dSJacob Faibussowitsch PetscCall(PetscRandomSetInterval(rand,-.001,.001)); 206*9566063dSJacob Faibussowitsch PetscCall(PetscRandomSetFromOptions(rand)); 207c4762a1bSJed Brown 208*9566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(da,0,&mx,&my,&mz,0,0,0,&dof,0,0,0,0,0)); 209c4762a1bSJed Brown /* For simplicity, this example only works on mx=my=mz */ 2102c71b3e2SJacob Faibussowitsch PetscCheckFalse(mx != my || mx != mz,PETSC_COMM_SELF,PETSC_ERR_SUP,"This example only works with mx %" PetscInt_FMT " = my %" PetscInt_FMT " = mz %" PetscInt_FMT,mx,my,mz); 211c4762a1bSJed Brown 212c4762a1bSJed Brown Hx = 1.0 / (PetscReal)(mx-1); Hy = 1.0 / (PetscReal)(my-1); Hz = 1.0 / (PetscReal)(mz-1); 213c4762a1bSJed Brown HxHydHz = Hx*Hy/Hz; HxHzdHy = Hx*Hz/Hy; HyHzdHx = Hy*Hz/Hx; 214c4762a1bSJed Brown 215*9566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(2*dof*dof+1,&v)); 216c4762a1bSJed Brown v_neighbor = v + dof*dof; 217*9566063dSJacob Faibussowitsch PetscCall(PetscArrayzero(v,2*dof*dof+1)); 218c4762a1bSJed Brown k3 = 0; 219c4762a1bSJed Brown for (k1=0; k1<dof; k1++) { 220c4762a1bSJed Brown for (k2=0; k2<dof; k2++) { 221c4762a1bSJed Brown if (k1 == k2) { 222c4762a1bSJed Brown v[k3] = 2.0*(HxHydHz + HxHzdHy + HyHzdHx); 223c4762a1bSJed Brown v_neighbor[k3] = -HxHydHz; 224c4762a1bSJed Brown } else { 225*9566063dSJacob Faibussowitsch PetscCall(PetscRandomGetValue(rand,&r1)); 226*9566063dSJacob Faibussowitsch PetscCall(PetscRandomGetValue(rand,&r2)); 227c4762a1bSJed Brown 228c4762a1bSJed Brown v[k3] = r1; 229c4762a1bSJed Brown v_neighbor[k3] = r2; 230c4762a1bSJed Brown } 231c4762a1bSJed Brown k3++; 232c4762a1bSJed Brown } 233c4762a1bSJed Brown } 234*9566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(da,&xs,&ys,&zs,&xm,&ym,&zm)); 235c4762a1bSJed Brown 236c4762a1bSJed Brown for (k=zs; k<zs+zm; k++) { 237c4762a1bSJed Brown for (j=ys; j<ys+ym; j++) { 238c4762a1bSJed Brown for (i=xs; i<xs+xm; i++) { 239c4762a1bSJed Brown row.i = i; row.j = j; row.k = k; 240a5b23f4aSJose E. Roman if (i==0 || j==0 || k==0 || i==mx-1 || j==my-1 || k==mz-1) { /* boundary points */ 241*9566063dSJacob Faibussowitsch PetscCall(MatSetValuesBlockedStencil(B,1,&row,1,&row,v,INSERT_VALUES)); 242c4762a1bSJed Brown } else { /* interior points */ 243c4762a1bSJed Brown /* center */ 244c4762a1bSJed Brown col.i = i; col.j = j; col.k = k; 245*9566063dSJacob Faibussowitsch PetscCall(MatSetValuesBlockedStencil(B,1,&row,1,&col,v,INSERT_VALUES)); 246c4762a1bSJed Brown 247c4762a1bSJed Brown /* x neighbors */ 248c4762a1bSJed Brown col.i = i-1; col.j = j; col.k = k; 249*9566063dSJacob Faibussowitsch PetscCall(MatSetValuesBlockedStencil(B,1,&row,1,&col,v_neighbor,INSERT_VALUES)); 250c4762a1bSJed Brown col.i = i+1; col.j = j; col.k = k; 251*9566063dSJacob Faibussowitsch PetscCall(MatSetValuesBlockedStencil(B,1,&row,1,&col,v_neighbor,INSERT_VALUES)); 252c4762a1bSJed Brown 253c4762a1bSJed Brown /* y neighbors */ 254c4762a1bSJed Brown col.i = i; col.j = j-1; col.k = k; 255*9566063dSJacob Faibussowitsch PetscCall(MatSetValuesBlockedStencil(B,1,&row,1,&col,v_neighbor,INSERT_VALUES)); 256c4762a1bSJed Brown col.i = i; col.j = j+1; col.k = k; 257*9566063dSJacob Faibussowitsch PetscCall(MatSetValuesBlockedStencil(B,1,&row,1,&col,v_neighbor,INSERT_VALUES)); 258c4762a1bSJed Brown 259c4762a1bSJed Brown /* z neighbors */ 260c4762a1bSJed Brown col.i = i; col.j = j; col.k = k-1; 261*9566063dSJacob Faibussowitsch PetscCall(MatSetValuesBlockedStencil(B,1,&row,1,&col,v_neighbor,INSERT_VALUES)); 262c4762a1bSJed Brown col.i = i; col.j = j; col.k = k+1; 263*9566063dSJacob Faibussowitsch PetscCall(MatSetValuesBlockedStencil(B,1,&row,1,&col,v_neighbor,INSERT_VALUES)); 264c4762a1bSJed Brown } 265c4762a1bSJed Brown } 266c4762a1bSJed Brown } 267c4762a1bSJed Brown } 268*9566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY)); 269*9566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY)); 270*9566063dSJacob Faibussowitsch PetscCall(PetscFree(v)); 271*9566063dSJacob Faibussowitsch PetscCall(PetscRandomDestroy(&rand)); 272c4762a1bSJed Brown PetscFunctionReturn(0); 273c4762a1bSJed Brown } 274c4762a1bSJed Brown 275c4762a1bSJed Brown /*TEST 276c4762a1bSJed Brown 277c4762a1bSJed Brown test: 278c4762a1bSJed Brown args: -dm_mat_type aij -dof 1 279c4762a1bSJed Brown output_file: output/ex129.out 280c4762a1bSJed Brown 281c4762a1bSJed Brown test: 282c4762a1bSJed Brown suffix: 2 283c4762a1bSJed Brown args: -dm_mat_type aij -dof 1 -inplacelu 284c4762a1bSJed Brown output_file: output/ex129.out 285c4762a1bSJed Brown 286c4762a1bSJed Brown TEST*/ 287