/* Laplacian in 3D. Use for testing MatSolve routines. Modeled by the partial differential equation - Laplacian u = 1,0 < x,y,z < 1, with boundary conditions u = 1 for x = 0, x = 1, y = 0, y = 1, z = 0, z = 1. */ static char help[] = "This example is for testing different MatSolve routines :MatSolve(), MatSolveAdd(), MatSolveTranspose(), MatSolveTransposeAdd(), and MatMatSolve().\n\ Example usage: ./ex129 -mat_type aij -dof 2\n\n"; #include #include extern PetscErrorCode ComputeMatrix(DM,Mat); extern PetscErrorCode ComputeRHS(DM,Vec); extern PetscErrorCode ComputeRHSMatrix(PetscInt,PetscInt,Mat*); int main(int argc,char **args) { PetscErrorCode ierr; PetscMPIInt size; Vec x,b,y,b1; DM da; Mat A,F,RHS,X,C1; MatFactorInfo info; IS perm,iperm; PetscInt dof =1,M=8,m,n,nrhs; PetscScalar one = 1.0; PetscReal norm,tol = 1000*PETSC_MACHINE_EPSILON; PetscBool InplaceLU=PETSC_FALSE; ierr = PetscInitialize(&argc,&args,(char*)0,help);if (ierr) return ierr; ierr = MPI_Comm_size(PETSC_COMM_WORLD,&size);CHKERRMPI(ierr); if (size != 1) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_WRONG_MPI_SIZE,"This is a uniprocessor example only\n"); ierr = PetscOptionsGetInt(NULL,NULL,"-dof",&dof,NULL);CHKERRQ(ierr); ierr = PetscOptionsGetInt(NULL,NULL,"-M",&M,NULL);CHKERRQ(ierr); ierr = DMDACreate(PETSC_COMM_WORLD,&da);CHKERRQ(ierr); ierr = DMSetDimension(da,3);CHKERRQ(ierr); ierr = DMDASetBoundaryType(da,DM_BOUNDARY_NONE,DM_BOUNDARY_NONE,DM_BOUNDARY_NONE);CHKERRQ(ierr); ierr = DMDASetStencilType(da,DMDA_STENCIL_STAR);CHKERRQ(ierr); ierr = DMDASetSizes(da,M,M,M);CHKERRQ(ierr); ierr = DMDASetNumProcs(da,PETSC_DECIDE,PETSC_DECIDE,PETSC_DECIDE);CHKERRQ(ierr); ierr = DMDASetDof(da,dof);CHKERRQ(ierr); ierr = DMDASetStencilWidth(da,1);CHKERRQ(ierr); ierr = DMDASetOwnershipRanges(da,NULL,NULL,NULL);CHKERRQ(ierr); ierr = DMSetMatType(da,MATBAIJ);CHKERRQ(ierr); ierr = DMSetFromOptions(da);CHKERRQ(ierr); ierr = DMSetUp(da);CHKERRQ(ierr); ierr = DMCreateGlobalVector(da,&x);CHKERRQ(ierr); ierr = DMCreateGlobalVector(da,&b);CHKERRQ(ierr); ierr = VecDuplicate(b,&y);CHKERRQ(ierr); ierr = ComputeRHS(da,b);CHKERRQ(ierr); ierr = VecSet(y,one);CHKERRQ(ierr); ierr = DMCreateMatrix(da,&A);CHKERRQ(ierr); ierr = ComputeMatrix(da,A);CHKERRQ(ierr); ierr = MatGetSize(A,&m,&n);CHKERRQ(ierr); nrhs = 2; ierr = PetscOptionsGetInt(NULL,NULL,"-nrhs",&nrhs,NULL);CHKERRQ(ierr); ierr = ComputeRHSMatrix(m,nrhs,&RHS);CHKERRQ(ierr); ierr = MatDuplicate(RHS,MAT_DO_NOT_COPY_VALUES,&X);CHKERRQ(ierr); ierr = MatGetOrdering(A,MATORDERINGND,&perm,&iperm);CHKERRQ(ierr); ierr = PetscOptionsGetBool(NULL,NULL,"-inplacelu",&InplaceLU,NULL);CHKERRQ(ierr); ierr = MatFactorInfoInitialize(&info);CHKERRQ(ierr); if (!InplaceLU) { ierr = MatGetFactor(A,MATSOLVERPETSC,MAT_FACTOR_LU,&F);CHKERRQ(ierr); info.fill = 5.0; ierr = MatLUFactorSymbolic(F,A,perm,iperm,&info);CHKERRQ(ierr); ierr = MatLUFactorNumeric(F,A,&info);CHKERRQ(ierr); } else { /* Test inplace factorization */ ierr = MatDuplicate(A,MAT_COPY_VALUES,&F);CHKERRQ(ierr); ierr = MatLUFactor(F,perm,iperm,&info);CHKERRQ(ierr); } ierr = VecDuplicate(y,&b1);CHKERRQ(ierr); /* MatSolve */ ierr = MatSolve(F,b,x);CHKERRQ(ierr); ierr = MatMult(A,x,b1);CHKERRQ(ierr); ierr = VecAXPY(b1,-1.0,b);CHKERRQ(ierr); ierr = VecNorm(b1,NORM_2,&norm);CHKERRQ(ierr); if (norm > tol) { ierr = PetscPrintf(PETSC_COMM_WORLD,"MatSolve : Error of norm %g\n",(double)norm);CHKERRQ(ierr); } /* MatSolveTranspose */ ierr = MatSolveTranspose(F,b,x);CHKERRQ(ierr); ierr = MatMultTranspose(A,x,b1);CHKERRQ(ierr); ierr = VecAXPY(b1,-1.0,b);CHKERRQ(ierr); ierr = VecNorm(b1,NORM_2,&norm);CHKERRQ(ierr); if (norm > tol) { ierr = PetscPrintf(PETSC_COMM_WORLD,"MatSolveTranspose : Error of norm %g\n",(double)norm);CHKERRQ(ierr); } /* MatSolveAdd */ ierr = MatSolveAdd(F,b,y,x);CHKERRQ(ierr); ierr = MatMult(A,y,b1);CHKERRQ(ierr); ierr = VecScale(b1,-1.0);CHKERRQ(ierr); ierr = MatMultAdd(A,x,b1,b1);CHKERRQ(ierr); ierr = VecAXPY(b1,-1.0,b);CHKERRQ(ierr); ierr = VecNorm(b1,NORM_2,&norm);CHKERRQ(ierr); if (norm > tol) { ierr = PetscPrintf(PETSC_COMM_WORLD,"MatSolveAdd : Error of norm %g\n",(double)norm);CHKERRQ(ierr); } /* MatSolveTransposeAdd */ ierr = MatSolveTransposeAdd(F,b,y,x);CHKERRQ(ierr); ierr = MatMultTranspose(A,y,b1);CHKERRQ(ierr); ierr = VecScale(b1,-1.0);CHKERRQ(ierr); ierr = MatMultTransposeAdd(A,x,b1,b1);CHKERRQ(ierr); ierr = VecAXPY(b1,-1.0,b);CHKERRQ(ierr); ierr = VecNorm(b1,NORM_2,&norm);CHKERRQ(ierr); if (norm > tol) { ierr = PetscPrintf(PETSC_COMM_WORLD,"MatSolveTransposeAdd : Error of norm %g\n",(double)norm);CHKERRQ(ierr); } /* MatMatSolve */ ierr = MatMatSolve(F,RHS,X);CHKERRQ(ierr); ierr = MatMatMult(A,X,MAT_INITIAL_MATRIX,2.0,&C1);CHKERRQ(ierr); ierr = MatAXPY(C1,-1.0,RHS,SAME_NONZERO_PATTERN);CHKERRQ(ierr); ierr = MatNorm(C1,NORM_FROBENIUS,&norm);CHKERRQ(ierr); if (norm > tol) { ierr = PetscPrintf(PETSC_COMM_WORLD,"MatMatSolve : Error of norm %g\n",(double)norm);CHKERRQ(ierr); } ierr = VecDestroy(&x);CHKERRQ(ierr); ierr = VecDestroy(&b);CHKERRQ(ierr); ierr = VecDestroy(&b1);CHKERRQ(ierr); ierr = VecDestroy(&y);CHKERRQ(ierr); ierr = MatDestroy(&A);CHKERRQ(ierr); ierr = MatDestroy(&F);CHKERRQ(ierr); ierr = MatDestroy(&RHS);CHKERRQ(ierr); ierr = MatDestroy(&C1);CHKERRQ(ierr); ierr = MatDestroy(&X);CHKERRQ(ierr); ierr = ISDestroy(&perm);CHKERRQ(ierr); ierr = ISDestroy(&iperm);CHKERRQ(ierr); ierr = DMDestroy(&da);CHKERRQ(ierr); ierr = PetscFinalize(); return ierr; } PetscErrorCode ComputeRHS(DM da,Vec b) { PetscErrorCode ierr; PetscInt mx,my,mz; PetscScalar h; PetscFunctionBegin; ierr = DMDAGetInfo(da,0,&mx,&my,&mz,0,0,0,0,0,0,0,0,0);CHKERRQ(ierr); h = 1.0/((mx-1)*(my-1)*(mz-1)); ierr = VecSet(b,h);CHKERRQ(ierr); PetscFunctionReturn(0); } PetscErrorCode ComputeRHSMatrix(PetscInt m,PetscInt nrhs,Mat *C) { PetscErrorCode ierr; PetscRandom rand; Mat RHS; PetscScalar *array,rval; PetscInt i,k; PetscFunctionBegin; ierr = MatCreate(PETSC_COMM_WORLD,&RHS);CHKERRQ(ierr); ierr = MatSetSizes(RHS,m,PETSC_DECIDE,PETSC_DECIDE,nrhs);CHKERRQ(ierr); ierr = MatSetType(RHS,MATSEQDENSE);CHKERRQ(ierr); ierr = MatSetUp(RHS);CHKERRQ(ierr); ierr = PetscRandomCreate(PETSC_COMM_WORLD,&rand);CHKERRQ(ierr); ierr = PetscRandomSetFromOptions(rand);CHKERRQ(ierr); ierr = MatDenseGetArray(RHS,&array);CHKERRQ(ierr); for (i=0; i 1) { for (k=1; k