1be1d678aSKris Buschelman 2b4319ba4SBarry Smith /* 3b4319ba4SBarry Smith Creates a matrix class for using the Neumann-Neumann type preconditioners. 4b4319ba4SBarry Smith This stores the matrices in globally unassembled form. Each processor 5b4319ba4SBarry Smith assembles only its local Neumann problem and the parallel matrix vector 6b4319ba4SBarry Smith product is handled "implicitly". 7b4319ba4SBarry Smith 8b4319ba4SBarry Smith We provide: 9b4319ba4SBarry Smith MatMult() 10b4319ba4SBarry Smith 11b4319ba4SBarry Smith Currently this allows for only one subdomain per processor. 12b4319ba4SBarry Smith 13b4319ba4SBarry Smith */ 14b4319ba4SBarry Smith 15c6db04a5SJed Brown #include <../src/mat/impls/is/matis.h> /*I "petscmat.h" I*/ 16b4319ba4SBarry Smith 17b4319ba4SBarry Smith #undef __FUNCT__ 18b4319ba4SBarry Smith #define __FUNCT__ "MatDestroy_IS" 19dfbe8321SBarry Smith PetscErrorCode MatDestroy_IS(Mat A) 20b4319ba4SBarry Smith { 21dfbe8321SBarry Smith PetscErrorCode ierr; 22b4319ba4SBarry Smith Mat_IS *b = (Mat_IS*)A->data; 23b4319ba4SBarry Smith 24b4319ba4SBarry Smith PetscFunctionBegin; 256bf464f9SBarry Smith ierr = MatDestroy(&b->A);CHKERRQ(ierr); 266bf464f9SBarry Smith ierr = VecScatterDestroy(&b->ctx);CHKERRQ(ierr); 276bf464f9SBarry Smith ierr = VecDestroy(&b->x);CHKERRQ(ierr); 286bf464f9SBarry Smith ierr = VecDestroy(&b->y);CHKERRQ(ierr); 296bf464f9SBarry Smith ierr = ISLocalToGlobalMappingDestroy(&b->mapping);CHKERRQ(ierr); 30bf0cc555SLisandro Dalcin ierr = PetscFree(A->data);CHKERRQ(ierr); 31dbd8c25aSHong Zhang ierr = PetscObjectChangeTypeName((PetscObject)A,0);CHKERRQ(ierr); 320298fd71SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)A,"MatISGetLocalMat_C","",NULL);CHKERRQ(ierr); 33b4319ba4SBarry Smith PetscFunctionReturn(0); 34b4319ba4SBarry Smith } 35b4319ba4SBarry Smith 36b4319ba4SBarry Smith #undef __FUNCT__ 37b4319ba4SBarry Smith #define __FUNCT__ "MatMult_IS" 38dfbe8321SBarry Smith PetscErrorCode MatMult_IS(Mat A,Vec x,Vec y) 39b4319ba4SBarry Smith { 40dfbe8321SBarry Smith PetscErrorCode ierr; 41b4319ba4SBarry Smith Mat_IS *is = (Mat_IS*)A->data; 42b4319ba4SBarry Smith PetscScalar zero = 0.0; 43b4319ba4SBarry Smith 44b4319ba4SBarry Smith PetscFunctionBegin; 45b4319ba4SBarry Smith /* scatter the global vector x into the local work vector */ 46ca9f406cSSatish Balay ierr = VecScatterBegin(is->ctx,x,is->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 47ca9f406cSSatish Balay ierr = VecScatterEnd(is->ctx,x,is->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 48b4319ba4SBarry Smith 49b4319ba4SBarry Smith /* multiply the local matrix */ 50b4319ba4SBarry Smith ierr = MatMult(is->A,is->x,is->y);CHKERRQ(ierr); 51b4319ba4SBarry Smith 52b4319ba4SBarry Smith /* scatter product back into global memory */ 532dcb1b2aSMatthew Knepley ierr = VecSet(y,zero);CHKERRQ(ierr); 54ca9f406cSSatish Balay ierr = VecScatterBegin(is->ctx,is->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 55ca9f406cSSatish Balay ierr = VecScatterEnd(is->ctx,is->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 56b4319ba4SBarry Smith PetscFunctionReturn(0); 57b4319ba4SBarry Smith } 58b4319ba4SBarry Smith 59b4319ba4SBarry Smith #undef __FUNCT__ 602e74eeadSLisandro Dalcin #define __FUNCT__ "MatMultAdd_IS" 612e74eeadSLisandro Dalcin static PetscErrorCode MatMultAdd_IS(Mat A,Vec v1,Vec v2,Vec v3) 622e74eeadSLisandro Dalcin { 63650997f4SStefano Zampini Vec temp_vec; 642e74eeadSLisandro Dalcin PetscErrorCode ierr; 652e74eeadSLisandro Dalcin 662e74eeadSLisandro Dalcin PetscFunctionBegin; /* v3 = v2 + A * v1.*/ 67650997f4SStefano Zampini if (v3 != v2) { 68650997f4SStefano Zampini ierr = MatMult(A,v1,v3);CHKERRQ(ierr); 69650997f4SStefano Zampini ierr = VecAXPY(v3,1.0,v2);CHKERRQ(ierr); 70650997f4SStefano Zampini } else { 71650997f4SStefano Zampini ierr = VecDuplicate(v2,&temp_vec);CHKERRQ(ierr); 72650997f4SStefano Zampini ierr = MatMult(A,v1,temp_vec);CHKERRQ(ierr); 73650997f4SStefano Zampini ierr = VecAXPY(temp_vec,1.0,v2);CHKERRQ(ierr); 74650997f4SStefano Zampini ierr = VecCopy(temp_vec,v3);CHKERRQ(ierr); 75650997f4SStefano Zampini ierr = VecDestroy(&temp_vec);CHKERRQ(ierr); 76650997f4SStefano Zampini } 772e74eeadSLisandro Dalcin PetscFunctionReturn(0); 782e74eeadSLisandro Dalcin } 792e74eeadSLisandro Dalcin 802e74eeadSLisandro Dalcin #undef __FUNCT__ 812e74eeadSLisandro Dalcin #define __FUNCT__ "MatMultTranspose_IS" 822e74eeadSLisandro Dalcin PetscErrorCode MatMultTranspose_IS(Mat A,Vec x,Vec y) 832e74eeadSLisandro Dalcin { 842e74eeadSLisandro Dalcin Mat_IS *is = (Mat_IS*)A->data; 852e74eeadSLisandro Dalcin PetscErrorCode ierr; 862e74eeadSLisandro Dalcin 872e74eeadSLisandro Dalcin PetscFunctionBegin; /* y = A' * x */ 882e74eeadSLisandro Dalcin /* scatter the global vector x into the local work vector */ 89ca9f406cSSatish Balay ierr = VecScatterBegin(is->ctx,x,is->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 90ca9f406cSSatish Balay ierr = VecScatterEnd(is->ctx,x,is->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 912e74eeadSLisandro Dalcin 922e74eeadSLisandro Dalcin /* multiply the local matrix */ 932e74eeadSLisandro Dalcin ierr = MatMultTranspose(is->A,is->x,is->y);CHKERRQ(ierr); 942e74eeadSLisandro Dalcin 952e74eeadSLisandro Dalcin /* scatter product back into global vector */ 962e74eeadSLisandro Dalcin ierr = VecSet(y,0);CHKERRQ(ierr); 97ca9f406cSSatish Balay ierr = VecScatterBegin(is->ctx,is->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 98ca9f406cSSatish Balay ierr = VecScatterEnd(is->ctx,is->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 992e74eeadSLisandro Dalcin PetscFunctionReturn(0); 1002e74eeadSLisandro Dalcin } 1012e74eeadSLisandro Dalcin 1022e74eeadSLisandro Dalcin #undef __FUNCT__ 1032e74eeadSLisandro Dalcin #define __FUNCT__ "MatMultTransposeAdd_IS" 1042e74eeadSLisandro Dalcin PetscErrorCode MatMultTransposeAdd_IS(Mat A,Vec v1,Vec v2,Vec v3) 1052e74eeadSLisandro Dalcin { 106650997f4SStefano Zampini Vec temp_vec; 1072e74eeadSLisandro Dalcin PetscErrorCode ierr; 1082e74eeadSLisandro Dalcin 1092e74eeadSLisandro Dalcin PetscFunctionBegin; /* v3 = v2 + A' * v1.*/ 110650997f4SStefano Zampini if (v3 != v2) { 111650997f4SStefano Zampini ierr = MatMultTranspose(A,v1,v3);CHKERRQ(ierr); 112650997f4SStefano Zampini ierr = VecAXPY(v3,1.0,v2);CHKERRQ(ierr); 113650997f4SStefano Zampini } else { 114650997f4SStefano Zampini ierr = VecDuplicate(v2,&temp_vec);CHKERRQ(ierr); 115650997f4SStefano Zampini ierr = MatMultTranspose(A,v1,temp_vec);CHKERRQ(ierr); 116650997f4SStefano Zampini ierr = VecAXPY(temp_vec,1.0,v2);CHKERRQ(ierr); 117650997f4SStefano Zampini ierr = VecCopy(temp_vec,v3);CHKERRQ(ierr); 118650997f4SStefano Zampini ierr = VecDestroy(&temp_vec);CHKERRQ(ierr); 119650997f4SStefano Zampini } 1202e74eeadSLisandro Dalcin PetscFunctionReturn(0); 1212e74eeadSLisandro Dalcin } 1222e74eeadSLisandro Dalcin 1232e74eeadSLisandro Dalcin #undef __FUNCT__ 124b4319ba4SBarry Smith #define __FUNCT__ "MatView_IS" 125dfbe8321SBarry Smith PetscErrorCode MatView_IS(Mat A,PetscViewer viewer) 126b4319ba4SBarry Smith { 127b4319ba4SBarry Smith Mat_IS *a = (Mat_IS*)A->data; 128dfbe8321SBarry Smith PetscErrorCode ierr; 129b4319ba4SBarry Smith PetscViewer sviewer; 130b4319ba4SBarry Smith 131b4319ba4SBarry Smith PetscFunctionBegin; 132b4319ba4SBarry Smith ierr = PetscViewerGetSingleton(viewer,&sviewer);CHKERRQ(ierr); 133b4319ba4SBarry Smith ierr = MatView(a->A,sviewer);CHKERRQ(ierr); 134b4319ba4SBarry Smith ierr = PetscViewerRestoreSingleton(viewer,&sviewer);CHKERRQ(ierr); 135b4319ba4SBarry Smith PetscFunctionReturn(0); 136b4319ba4SBarry Smith } 137b4319ba4SBarry Smith 138b4319ba4SBarry Smith #undef __FUNCT__ 139b4319ba4SBarry Smith #define __FUNCT__ "MatSetLocalToGlobalMapping_IS" 140784ac674SJed Brown PetscErrorCode MatSetLocalToGlobalMapping_IS(Mat A,ISLocalToGlobalMapping rmapping,ISLocalToGlobalMapping cmapping) 141b4319ba4SBarry Smith { 142dfbe8321SBarry Smith PetscErrorCode ierr; 1434e4c7dbeSStefano Zampini PetscInt n,bs; 144b4319ba4SBarry Smith Mat_IS *is = (Mat_IS*)A->data; 145b4319ba4SBarry Smith IS from,to; 146b4319ba4SBarry Smith Vec global; 147b4319ba4SBarry Smith 148b4319ba4SBarry Smith PetscFunctionBegin; 149e7e72b3dSBarry Smith if (is->mapping) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Mapping already set for matrix"); 150784ac674SJed Brown PetscCheckSameComm(A,1,rmapping,2); 151ce94432eSBarry Smith if (rmapping != cmapping) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_INCOMP,"MATIS requires the row and column mappings to be identical"); 152784ac674SJed Brown ierr = PetscObjectReference((PetscObject)rmapping);CHKERRQ(ierr); 1536bf464f9SBarry Smith ierr = ISLocalToGlobalMappingDestroy(&is->mapping);CHKERRQ(ierr); 154784ac674SJed Brown is->mapping = rmapping; 155b4319ba4SBarry Smith 156b4319ba4SBarry Smith /* Create the local matrix A */ 157784ac674SJed Brown ierr = ISLocalToGlobalMappingGetSize(rmapping,&n);CHKERRQ(ierr); 1584e4c7dbeSStefano Zampini ierr = MatGetBlockSize(A,&bs);CHKERRQ(ierr); 159f69a0ea3SMatthew Knepley ierr = MatCreate(PETSC_COMM_SELF,&is->A);CHKERRQ(ierr); 160f69a0ea3SMatthew Knepley ierr = MatSetSizes(is->A,n,n,n,n);CHKERRQ(ierr); 1614e4c7dbeSStefano Zampini ierr = MatSetBlockSize(is->A,bs);CHKERRQ(ierr); 162ff130e51SJed Brown ierr = MatSetOptionsPrefix(is->A,((PetscObject)A)->prefix);CHKERRQ(ierr); 163ff130e51SJed Brown ierr = MatAppendOptionsPrefix(is->A,"is_");CHKERRQ(ierr); 164b4319ba4SBarry Smith ierr = MatSetFromOptions(is->A);CHKERRQ(ierr); 165b4319ba4SBarry Smith 166b4319ba4SBarry Smith /* Create the local work vectors */ 1674e4c7dbeSStefano Zampini ierr = VecCreate(PETSC_COMM_SELF,&is->x);CHKERRQ(ierr); 1684e4c7dbeSStefano Zampini ierr = VecSetBlockSize(is->x,bs);CHKERRQ(ierr); 1694e4c7dbeSStefano Zampini ierr = VecSetSizes(is->x,n,n);CHKERRQ(ierr); 170ff130e51SJed Brown ierr = VecSetOptionsPrefix(is->x,((PetscObject)A)->prefix);CHKERRQ(ierr); 171ff130e51SJed Brown ierr = VecAppendOptionsPrefix(is->x,"is_");CHKERRQ(ierr); 1724e4c7dbeSStefano Zampini ierr = VecSetFromOptions(is->x);CHKERRQ(ierr); 173b4319ba4SBarry Smith ierr = VecDuplicate(is->x,&is->y);CHKERRQ(ierr); 174b4319ba4SBarry Smith 175b4319ba4SBarry Smith /* setup the global to local scatter */ 176b4319ba4SBarry Smith ierr = ISCreateStride(PETSC_COMM_SELF,n,0,1,&to);CHKERRQ(ierr); 177784ac674SJed Brown ierr = ISLocalToGlobalMappingApplyIS(rmapping,to,&from);CHKERRQ(ierr); 1780298fd71SBarry Smith ierr = MatGetVecs(A,&global,NULL);CHKERRQ(ierr); 179b4319ba4SBarry Smith ierr = VecScatterCreate(global,from,is->x,to,&is->ctx);CHKERRQ(ierr); 1806bf464f9SBarry Smith ierr = VecDestroy(&global);CHKERRQ(ierr); 1816bf464f9SBarry Smith ierr = ISDestroy(&to);CHKERRQ(ierr); 1826bf464f9SBarry Smith ierr = ISDestroy(&from);CHKERRQ(ierr); 183b4319ba4SBarry Smith PetscFunctionReturn(0); 184b4319ba4SBarry Smith } 185b4319ba4SBarry Smith 1862e74eeadSLisandro Dalcin #undef __FUNCT__ 1872e74eeadSLisandro Dalcin #define __FUNCT__ "MatSetValues_IS" 1882e74eeadSLisandro Dalcin PetscErrorCode MatSetValues_IS(Mat mat, PetscInt m,const PetscInt *rows, PetscInt n,const PetscInt *cols, const PetscScalar *values, InsertMode addv) 1892e74eeadSLisandro Dalcin { 1902e74eeadSLisandro Dalcin Mat_IS *is = (Mat_IS*)mat->data; 1912e74eeadSLisandro Dalcin PetscInt rows_l[2048],cols_l[2048]; 1922e74eeadSLisandro Dalcin PetscErrorCode ierr; 1932e74eeadSLisandro Dalcin 1942e74eeadSLisandro Dalcin PetscFunctionBegin; 1952e74eeadSLisandro Dalcin #if defined(PETSC_USE_DEBUG) 196f23aa3ddSBarry Smith if (m > 2048 || n > 2048) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_SUP,"Number of row/column indices must be <= 2048: they are %D %D",m,n); 1972e74eeadSLisandro Dalcin #endif 1982e74eeadSLisandro Dalcin ierr = ISG2LMapApply(is->mapping,m,rows,rows_l);CHKERRQ(ierr); 1992e74eeadSLisandro Dalcin ierr = ISG2LMapApply(is->mapping,n,cols,cols_l);CHKERRQ(ierr); 2002e74eeadSLisandro Dalcin ierr = MatSetValues(is->A,m,rows_l,n,cols_l,values,addv);CHKERRQ(ierr); 2012e74eeadSLisandro Dalcin PetscFunctionReturn(0); 2022e74eeadSLisandro Dalcin } 2032e74eeadSLisandro Dalcin 2042e74eeadSLisandro Dalcin #undef ISG2LMapSetUp 2052e74eeadSLisandro Dalcin #undef ISG2LMapApply 206b4319ba4SBarry Smith 207b4319ba4SBarry Smith #undef __FUNCT__ 208b4319ba4SBarry Smith #define __FUNCT__ "MatSetValuesLocal_IS" 20913f74950SBarry Smith PetscErrorCode MatSetValuesLocal_IS(Mat A,PetscInt m,const PetscInt *rows, PetscInt n,const PetscInt *cols,const PetscScalar *values,InsertMode addv) 210b4319ba4SBarry Smith { 211dfbe8321SBarry Smith PetscErrorCode ierr; 212b4319ba4SBarry Smith Mat_IS *is = (Mat_IS*)A->data; 213b4319ba4SBarry Smith 214b4319ba4SBarry Smith PetscFunctionBegin; 215b4319ba4SBarry Smith ierr = MatSetValues(is->A,m,rows,n,cols,values,addv);CHKERRQ(ierr); 216b4319ba4SBarry Smith PetscFunctionReturn(0); 217b4319ba4SBarry Smith } 218b4319ba4SBarry Smith 219b4319ba4SBarry Smith #undef __FUNCT__ 2202e74eeadSLisandro Dalcin #define __FUNCT__ "MatZeroRows_IS" 2212b40b63fSBarry Smith PetscErrorCode MatZeroRows_IS(Mat A,PetscInt n,const PetscInt rows[],PetscScalar diag,Vec x,Vec b) 2222e74eeadSLisandro Dalcin { 2232e74eeadSLisandro Dalcin Mat_IS *is = (Mat_IS*)A->data; 2240298fd71SBarry Smith PetscInt n_l = 0, *rows_l = NULL; 2252e74eeadSLisandro Dalcin PetscErrorCode ierr; 2262e74eeadSLisandro Dalcin 2272e74eeadSLisandro Dalcin PetscFunctionBegin; 228ce94432eSBarry Smith if (x && b) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"No support"); 2292e74eeadSLisandro Dalcin if (n) { 2302e74eeadSLisandro Dalcin ierr = PetscMalloc(n*sizeof(PetscInt),&rows_l);CHKERRQ(ierr); 2312e74eeadSLisandro Dalcin ierr = ISGlobalToLocalMappingApply(is->mapping,IS_GTOLM_DROP,n,rows,&n_l,rows_l);CHKERRQ(ierr); 2322e74eeadSLisandro Dalcin } 2332b40b63fSBarry Smith ierr = MatZeroRowsLocal(A,n_l,rows_l,diag,x,b);CHKERRQ(ierr); 234c31cb41cSBarry Smith ierr = PetscFree(rows_l);CHKERRQ(ierr); 2352e74eeadSLisandro Dalcin PetscFunctionReturn(0); 2362e74eeadSLisandro Dalcin } 2372e74eeadSLisandro Dalcin 2382e74eeadSLisandro Dalcin #undef __FUNCT__ 239b4319ba4SBarry Smith #define __FUNCT__ "MatZeroRowsLocal_IS" 2402b40b63fSBarry Smith PetscErrorCode MatZeroRowsLocal_IS(Mat A,PetscInt n,const PetscInt rows[],PetscScalar diag,Vec x,Vec b) 241b4319ba4SBarry Smith { 242b4319ba4SBarry Smith Mat_IS *is = (Mat_IS*)A->data; 243dfbe8321SBarry Smith PetscErrorCode ierr; 244f4df32b1SMatthew Knepley PetscInt i; 245b4319ba4SBarry Smith PetscScalar *array; 246b4319ba4SBarry Smith 247b4319ba4SBarry Smith PetscFunctionBegin; 248ce94432eSBarry Smith if (x && b) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"No support"); 249b4319ba4SBarry Smith { 250b4319ba4SBarry Smith /* 251b4319ba4SBarry Smith Set up is->x as a "counting vector". This is in order to MatMult_IS 252b4319ba4SBarry Smith work properly in the interface nodes. 253b4319ba4SBarry Smith */ 254b4319ba4SBarry Smith Vec counter; 255b4319ba4SBarry Smith PetscScalar one=1.0, zero=0.0; 2560298fd71SBarry Smith ierr = MatGetVecs(A,&counter,NULL);CHKERRQ(ierr); 2572dcb1b2aSMatthew Knepley ierr = VecSet(counter,zero);CHKERRQ(ierr); 2582dcb1b2aSMatthew Knepley ierr = VecSet(is->x,one);CHKERRQ(ierr); 259ca9f406cSSatish Balay ierr = VecScatterBegin(is->ctx,is->x,counter,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 260ca9f406cSSatish Balay ierr = VecScatterEnd (is->ctx,is->x,counter,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 261ca9f406cSSatish Balay ierr = VecScatterBegin(is->ctx,counter,is->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 262ca9f406cSSatish Balay ierr = VecScatterEnd (is->ctx,counter,is->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 2636bf464f9SBarry Smith ierr = VecDestroy(&counter);CHKERRQ(ierr); 264b4319ba4SBarry Smith } 265958c9bccSBarry Smith if (!n) { 266b4319ba4SBarry Smith is->pure_neumann = PETSC_TRUE; 267b4319ba4SBarry Smith } else { 268b4319ba4SBarry Smith is->pure_neumann = PETSC_FALSE; 2692205254eSKarl Rupp 270b4319ba4SBarry Smith ierr = VecGetArray(is->x,&array);CHKERRQ(ierr); 2712b40b63fSBarry Smith ierr = MatZeroRows(is->A,n,rows,diag,0,0);CHKERRQ(ierr); 272b4319ba4SBarry Smith for (i=0; i<n; i++) { 273f4df32b1SMatthew Knepley ierr = MatSetValue(is->A,rows[i],rows[i],diag/(array[rows[i]]),INSERT_VALUES);CHKERRQ(ierr); 274b4319ba4SBarry Smith } 275b4319ba4SBarry Smith ierr = MatAssemblyBegin(is->A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 276b4319ba4SBarry Smith ierr = MatAssemblyEnd (is->A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 277b4319ba4SBarry Smith ierr = VecRestoreArray(is->x,&array);CHKERRQ(ierr); 278b4319ba4SBarry Smith } 279b4319ba4SBarry Smith PetscFunctionReturn(0); 280b4319ba4SBarry Smith } 281b4319ba4SBarry Smith 282b4319ba4SBarry Smith #undef __FUNCT__ 283b4319ba4SBarry Smith #define __FUNCT__ "MatAssemblyBegin_IS" 284dfbe8321SBarry Smith PetscErrorCode MatAssemblyBegin_IS(Mat A,MatAssemblyType type) 285b4319ba4SBarry Smith { 286b4319ba4SBarry Smith Mat_IS *is = (Mat_IS*)A->data; 287dfbe8321SBarry Smith PetscErrorCode ierr; 288b4319ba4SBarry Smith 289b4319ba4SBarry Smith PetscFunctionBegin; 290b4319ba4SBarry Smith ierr = MatAssemblyBegin(is->A,type);CHKERRQ(ierr); 291b4319ba4SBarry Smith PetscFunctionReturn(0); 292b4319ba4SBarry Smith } 293b4319ba4SBarry Smith 294b4319ba4SBarry Smith #undef __FUNCT__ 295b4319ba4SBarry Smith #define __FUNCT__ "MatAssemblyEnd_IS" 296dfbe8321SBarry Smith PetscErrorCode MatAssemblyEnd_IS(Mat A,MatAssemblyType type) 297b4319ba4SBarry Smith { 298b4319ba4SBarry Smith Mat_IS *is = (Mat_IS*)A->data; 299dfbe8321SBarry Smith PetscErrorCode ierr; 300b4319ba4SBarry Smith 301b4319ba4SBarry Smith PetscFunctionBegin; 302b4319ba4SBarry Smith ierr = MatAssemblyEnd(is->A,type);CHKERRQ(ierr); 303b4319ba4SBarry Smith PetscFunctionReturn(0); 304b4319ba4SBarry Smith } 305b4319ba4SBarry Smith 306b4319ba4SBarry Smith #undef __FUNCT__ 307b4319ba4SBarry Smith #define __FUNCT__ "MatISGetLocalMat_IS" 3087087cfbeSBarry Smith PetscErrorCode MatISGetLocalMat_IS(Mat mat,Mat *local) 309b4319ba4SBarry Smith { 310b4319ba4SBarry Smith Mat_IS *is = (Mat_IS*)mat->data; 311b4319ba4SBarry Smith 312b4319ba4SBarry Smith PetscFunctionBegin; 313b4319ba4SBarry Smith *local = is->A; 314b4319ba4SBarry Smith PetscFunctionReturn(0); 315b4319ba4SBarry Smith } 316b4319ba4SBarry Smith 317b4319ba4SBarry Smith #undef __FUNCT__ 318b4319ba4SBarry Smith #define __FUNCT__ "MatISGetLocalMat" 319b4319ba4SBarry Smith /*@ 320b4319ba4SBarry Smith MatISGetLocalMat - Gets the local matrix stored inside a MATIS matrix. 321b4319ba4SBarry Smith 322b4319ba4SBarry Smith Input Parameter: 323b4319ba4SBarry Smith . mat - the matrix 324b4319ba4SBarry Smith 325b4319ba4SBarry Smith Output Parameter: 326b4319ba4SBarry Smith . local - the local matrix usually MATSEQAIJ 327b4319ba4SBarry Smith 328b4319ba4SBarry Smith Level: advanced 329b4319ba4SBarry Smith 330b4319ba4SBarry Smith Notes: 331b4319ba4SBarry Smith This can be called if you have precomputed the nonzero structure of the 332b4319ba4SBarry Smith matrix and want to provide it to the inner matrix object to improve the performance 333b4319ba4SBarry Smith of the MatSetValues() operation. 334b4319ba4SBarry Smith 335b4319ba4SBarry Smith .seealso: MATIS 336b4319ba4SBarry Smith @*/ 3377087cfbeSBarry Smith PetscErrorCode MatISGetLocalMat(Mat mat,Mat *local) 338b4319ba4SBarry Smith { 3394ac538c5SBarry Smith PetscErrorCode ierr; 340b4319ba4SBarry Smith 341b4319ba4SBarry Smith PetscFunctionBegin; 3420700a824SBarry Smith PetscValidHeaderSpecific(mat,MAT_CLASSID,1); 343b4319ba4SBarry Smith PetscValidPointer(local,2); 3444ac538c5SBarry Smith ierr = PetscUseMethod(mat,"MatISGetLocalMat_C",(Mat,Mat*),(mat,local));CHKERRQ(ierr); 345b4319ba4SBarry Smith PetscFunctionReturn(0); 346b4319ba4SBarry Smith } 347b4319ba4SBarry Smith 3483b03a366Sstefano_zampini #undef __FUNCT__ 3493b03a366Sstefano_zampini #define __FUNCT__ "MatISSetLocalMat_IS" 3503b03a366Sstefano_zampini PetscErrorCode MatISSetLocalMat_IS(Mat mat,Mat local) 3513b03a366Sstefano_zampini { 3523b03a366Sstefano_zampini Mat_IS *is = (Mat_IS*)mat->data; 3533b03a366Sstefano_zampini PetscInt nrows,ncols,orows,ocols; 3543b03a366Sstefano_zampini PetscErrorCode ierr; 3553b03a366Sstefano_zampini 3563b03a366Sstefano_zampini PetscFunctionBegin; 3574e4c7dbeSStefano Zampini if (is->A) { 3583b03a366Sstefano_zampini ierr = MatGetSize(is->A,&orows,&ocols);CHKERRQ(ierr); 3593b03a366Sstefano_zampini ierr = MatGetSize(local,&nrows,&ncols);CHKERRQ(ierr); 360f23aa3ddSBarry Smith if (orows != nrows || ocols != ncols) SETERRQ4(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Local MATIS matrix should be of size %dx%d (you passed a %dx%d matrix)\n",orows,ocols,nrows,ncols); 3614e4c7dbeSStefano Zampini } 3623b03a366Sstefano_zampini ierr = PetscObjectReference((PetscObject)local);CHKERRQ(ierr); 3633b03a366Sstefano_zampini ierr = MatDestroy(&is->A);CHKERRQ(ierr); 3643b03a366Sstefano_zampini is->A = local; 3653b03a366Sstefano_zampini PetscFunctionReturn(0); 3663b03a366Sstefano_zampini } 3673b03a366Sstefano_zampini 3683b03a366Sstefano_zampini #undef __FUNCT__ 3693b03a366Sstefano_zampini #define __FUNCT__ "MatISSetLocalMat" 3703b03a366Sstefano_zampini /*@ 3713b03a366Sstefano_zampini MatISSetLocalMat - Set the local matrix stored inside a MATIS matrix. 3723b03a366Sstefano_zampini 3733b03a366Sstefano_zampini Input Parameter: 3743b03a366Sstefano_zampini . mat - the matrix 3753b03a366Sstefano_zampini . local - the local matrix usually MATSEQAIJ 3763b03a366Sstefano_zampini 3773b03a366Sstefano_zampini Output Parameter: 3783b03a366Sstefano_zampini 3793b03a366Sstefano_zampini Level: advanced 3803b03a366Sstefano_zampini 3813b03a366Sstefano_zampini Notes: 3823b03a366Sstefano_zampini This can be called if you have precomputed the local matrix and 3833b03a366Sstefano_zampini want to provide it to the matrix object MATIS. 3843b03a366Sstefano_zampini 3853b03a366Sstefano_zampini .seealso: MATIS 3863b03a366Sstefano_zampini @*/ 3873b03a366Sstefano_zampini PetscErrorCode MatISSetLocalMat(Mat mat,Mat local) 3883b03a366Sstefano_zampini { 3893b03a366Sstefano_zampini PetscErrorCode ierr; 3903b03a366Sstefano_zampini 3913b03a366Sstefano_zampini PetscFunctionBegin; 3923b03a366Sstefano_zampini PetscValidHeaderSpecific(mat,MAT_CLASSID,1); 3933b03a366Sstefano_zampini ierr = PetscUseMethod(mat,"MatISSetLocalMat_C",(Mat,Mat),(mat,local));CHKERRQ(ierr); 3943b03a366Sstefano_zampini PetscFunctionReturn(0); 3953b03a366Sstefano_zampini } 3963b03a366Sstefano_zampini 3976726f965SBarry Smith #undef __FUNCT__ 3986726f965SBarry Smith #define __FUNCT__ "MatZeroEntries_IS" 3996726f965SBarry Smith PetscErrorCode MatZeroEntries_IS(Mat A) 4006726f965SBarry Smith { 4016726f965SBarry Smith Mat_IS *a = (Mat_IS*)A->data; 4026726f965SBarry Smith PetscErrorCode ierr; 4036726f965SBarry Smith 4046726f965SBarry Smith PetscFunctionBegin; 4056726f965SBarry Smith ierr = MatZeroEntries(a->A);CHKERRQ(ierr); 4066726f965SBarry Smith PetscFunctionReturn(0); 4076726f965SBarry Smith } 4086726f965SBarry Smith 4096726f965SBarry Smith #undef __FUNCT__ 4102e74eeadSLisandro Dalcin #define __FUNCT__ "MatScale_IS" 4112e74eeadSLisandro Dalcin PetscErrorCode MatScale_IS(Mat A,PetscScalar a) 4122e74eeadSLisandro Dalcin { 4132e74eeadSLisandro Dalcin Mat_IS *is = (Mat_IS*)A->data; 4142e74eeadSLisandro Dalcin PetscErrorCode ierr; 4152e74eeadSLisandro Dalcin 4162e74eeadSLisandro Dalcin PetscFunctionBegin; 4172e74eeadSLisandro Dalcin ierr = MatScale(is->A,a);CHKERRQ(ierr); 4182e74eeadSLisandro Dalcin PetscFunctionReturn(0); 4192e74eeadSLisandro Dalcin } 4202e74eeadSLisandro Dalcin 4212e74eeadSLisandro Dalcin #undef __FUNCT__ 4222e74eeadSLisandro Dalcin #define __FUNCT__ "MatGetDiagonal_IS" 4232e74eeadSLisandro Dalcin PetscErrorCode MatGetDiagonal_IS(Mat A, Vec v) 4242e74eeadSLisandro Dalcin { 4252e74eeadSLisandro Dalcin Mat_IS *is = (Mat_IS*)A->data; 4262e74eeadSLisandro Dalcin PetscErrorCode ierr; 4272e74eeadSLisandro Dalcin 4282e74eeadSLisandro Dalcin PetscFunctionBegin; 4292e74eeadSLisandro Dalcin /* get diagonal of the local matrix */ 4302e74eeadSLisandro Dalcin ierr = MatGetDiagonal(is->A,is->x);CHKERRQ(ierr); 4312e74eeadSLisandro Dalcin 4322e74eeadSLisandro Dalcin /* scatter diagonal back into global vector */ 4332e74eeadSLisandro Dalcin ierr = VecSet(v,0);CHKERRQ(ierr); 434ca9f406cSSatish Balay ierr = VecScatterBegin(is->ctx,is->x,v,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 435ca9f406cSSatish Balay ierr = VecScatterEnd (is->ctx,is->x,v,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 4362e74eeadSLisandro Dalcin PetscFunctionReturn(0); 4372e74eeadSLisandro Dalcin } 4382e74eeadSLisandro Dalcin 4392e74eeadSLisandro Dalcin #undef __FUNCT__ 4406726f965SBarry Smith #define __FUNCT__ "MatSetOption_IS" 441ace3abfcSBarry Smith PetscErrorCode MatSetOption_IS(Mat A,MatOption op,PetscBool flg) 4426726f965SBarry Smith { 4436726f965SBarry Smith Mat_IS *a = (Mat_IS*)A->data; 4446726f965SBarry Smith PetscErrorCode ierr; 4456726f965SBarry Smith 4466726f965SBarry Smith PetscFunctionBegin; 4474e0d8c25SBarry Smith ierr = MatSetOption(a->A,op,flg);CHKERRQ(ierr); 4486726f965SBarry Smith PetscFunctionReturn(0); 4496726f965SBarry Smith } 4506726f965SBarry Smith 451284134d9SBarry Smith #undef __FUNCT__ 452284134d9SBarry Smith #define __FUNCT__ "MatCreateIS" 453284134d9SBarry Smith /*@ 454284134d9SBarry Smith MatCreateIS - Creates a "process" unassmembled matrix, it is assembled on each 455284134d9SBarry Smith process but not across processes. 456284134d9SBarry Smith 457284134d9SBarry Smith Input Parameters: 458284134d9SBarry Smith + comm - MPI communicator that will share the matrix 4594e4c7dbeSStefano Zampini . bs - local and global block size of the matrix 460284134d9SBarry Smith . m,n,M,N - local and/or global sizes of the the left and right vector used in matrix vector products 461284134d9SBarry Smith - map - mapping that defines the global number for each local number 462284134d9SBarry Smith 463284134d9SBarry Smith Output Parameter: 464284134d9SBarry Smith . A - the resulting matrix 465284134d9SBarry Smith 4668e6c10adSSatish Balay Level: advanced 4678e6c10adSSatish Balay 468284134d9SBarry Smith Notes: See MATIS for more details 4698cda6cd7SBarry Smith m and n are NOT related to the size of the map, they are the size of the part of the vector owned 4708cda6cd7SBarry Smith by that process. m + nghosts (or n + nghosts) is the length of map since map maps all local points 4718cda6cd7SBarry Smith plus the ghost points to global indices. 472284134d9SBarry Smith 473284134d9SBarry Smith .seealso: MATIS, MatSetLocalToGlobalMapping() 474284134d9SBarry Smith @*/ 4754e4c7dbeSStefano Zampini PetscErrorCode MatCreateIS(MPI_Comm comm,PetscInt bs,PetscInt m,PetscInt n,PetscInt M,PetscInt N,ISLocalToGlobalMapping map,Mat *A) 476284134d9SBarry Smith { 477284134d9SBarry Smith PetscErrorCode ierr; 478284134d9SBarry Smith 479284134d9SBarry Smith PetscFunctionBegin; 480284134d9SBarry Smith ierr = MatCreate(comm,A);CHKERRQ(ierr); 481d3ee2243SStefano Zampini ierr = MatSetBlockSize(*A,bs);CHKERRQ(ierr); 482284134d9SBarry Smith ierr = MatSetSizes(*A,m,n,M,N);CHKERRQ(ierr); 483284134d9SBarry Smith ierr = MatSetType(*A,MATIS);CHKERRQ(ierr); 4843b03a366Sstefano_zampini ierr = MatSetUp(*A);CHKERRQ(ierr); 485784ac674SJed Brown ierr = MatSetLocalToGlobalMapping(*A,map,map);CHKERRQ(ierr); 486284134d9SBarry Smith PetscFunctionReturn(0); 487284134d9SBarry Smith } 488284134d9SBarry Smith 489b4319ba4SBarry Smith /*MC 490b4319ba4SBarry Smith MATIS - MATIS = "is" - A matrix type to be used for using the Neumann-Neumann type preconditioners. 491b4319ba4SBarry Smith This stores the matrices in globally unassembled form. Each processor 492b4319ba4SBarry Smith assembles only its local Neumann problem and the parallel matrix vector 493b4319ba4SBarry Smith product is handled "implicitly". 494b4319ba4SBarry Smith 495b4319ba4SBarry Smith Operations Provided: 4966726f965SBarry Smith + MatMult() 4972e74eeadSLisandro Dalcin . MatMultAdd() 4982e74eeadSLisandro Dalcin . MatMultTranspose() 4992e74eeadSLisandro Dalcin . MatMultTransposeAdd() 5006726f965SBarry Smith . MatZeroEntries() 5016726f965SBarry Smith . MatSetOption() 5022e74eeadSLisandro Dalcin . MatZeroRows() 5036726f965SBarry Smith . MatZeroRowsLocal() 5042e74eeadSLisandro Dalcin . MatSetValues() 5056726f965SBarry Smith . MatSetValuesLocal() 5062e74eeadSLisandro Dalcin . MatScale() 5072e74eeadSLisandro Dalcin . MatGetDiagonal() 5086726f965SBarry Smith - MatSetLocalToGlobalMapping() 509b4319ba4SBarry Smith 510b4319ba4SBarry Smith Options Database Keys: 511b4319ba4SBarry Smith . -mat_type is - sets the matrix type to "is" during a call to MatSetFromOptions() 512b4319ba4SBarry Smith 513b4319ba4SBarry Smith Notes: Options prefix for the inner matrix are given by -is_mat_xxx 514b4319ba4SBarry Smith 515b4319ba4SBarry Smith You must call MatSetLocalToGlobalMapping() before using this matrix type. 516b4319ba4SBarry Smith 517b4319ba4SBarry Smith You can do matrix preallocation on the local matrix after you obtain it with 518b4319ba4SBarry Smith MatISGetLocalMat() 519b4319ba4SBarry Smith 520b4319ba4SBarry Smith Level: advanced 521b4319ba4SBarry Smith 522b4319ba4SBarry Smith .seealso: PC, MatISGetLocalMat(), MatSetLocalToGlobalMapping() 523b4319ba4SBarry Smith 524b4319ba4SBarry Smith M*/ 525b4319ba4SBarry Smith 526b4319ba4SBarry Smith #undef __FUNCT__ 527b4319ba4SBarry Smith #define __FUNCT__ "MatCreate_IS" 528*8cc058d9SJed Brown PETSC_EXTERN PetscErrorCode MatCreate_IS(Mat A) 529b4319ba4SBarry Smith { 530dfbe8321SBarry Smith PetscErrorCode ierr; 531b4319ba4SBarry Smith Mat_IS *b; 532b4319ba4SBarry Smith 533b4319ba4SBarry Smith PetscFunctionBegin; 53438f2d2fdSLisandro Dalcin ierr = PetscNewLog(A,Mat_IS,&b);CHKERRQ(ierr); 535b4319ba4SBarry Smith A->data = (void*)b; 536b4319ba4SBarry Smith ierr = PetscMemzero(A->ops,sizeof(struct _MatOps));CHKERRQ(ierr); 537b4319ba4SBarry Smith 538b4319ba4SBarry Smith A->ops->mult = MatMult_IS; 5392e74eeadSLisandro Dalcin A->ops->multadd = MatMultAdd_IS; 5402e74eeadSLisandro Dalcin A->ops->multtranspose = MatMultTranspose_IS; 5412e74eeadSLisandro Dalcin A->ops->multtransposeadd = MatMultTransposeAdd_IS; 542b4319ba4SBarry Smith A->ops->destroy = MatDestroy_IS; 543b4319ba4SBarry Smith A->ops->setlocaltoglobalmapping = MatSetLocalToGlobalMapping_IS; 5442e74eeadSLisandro Dalcin A->ops->setvalues = MatSetValues_IS; 545b4319ba4SBarry Smith A->ops->setvalueslocal = MatSetValuesLocal_IS; 5462e74eeadSLisandro Dalcin A->ops->zerorows = MatZeroRows_IS; 547b4319ba4SBarry Smith A->ops->zerorowslocal = MatZeroRowsLocal_IS; 548b4319ba4SBarry Smith A->ops->assemblybegin = MatAssemblyBegin_IS; 549b4319ba4SBarry Smith A->ops->assemblyend = MatAssemblyEnd_IS; 550b4319ba4SBarry Smith A->ops->view = MatView_IS; 5516726f965SBarry Smith A->ops->zeroentries = MatZeroEntries_IS; 5522e74eeadSLisandro Dalcin A->ops->scale = MatScale_IS; 5532e74eeadSLisandro Dalcin A->ops->getdiagonal = MatGetDiagonal_IS; 5546726f965SBarry Smith A->ops->setoption = MatSetOption_IS; 555b4319ba4SBarry Smith 55626283091SBarry Smith ierr = PetscLayoutSetUp(A->rmap);CHKERRQ(ierr); 55726283091SBarry Smith ierr = PetscLayoutSetUp(A->cmap);CHKERRQ(ierr); 558b4319ba4SBarry Smith 559b4319ba4SBarry Smith b->A = 0; 560b4319ba4SBarry Smith b->ctx = 0; 561b4319ba4SBarry Smith b->x = 0; 562b4319ba4SBarry Smith b->y = 0; 56300de8ff0SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)A,"MatISGetLocalMat_C","MatISGetLocalMat_IS",MatISGetLocalMat_IS);CHKERRQ(ierr); 56400de8ff0SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)A,"MatISSetLocalMat_C","MatISSetLocalMat_IS",MatISSetLocalMat_IS);CHKERRQ(ierr); 56517667f90SBarry Smith ierr = PetscObjectChangeTypeName((PetscObject)A,MATIS);CHKERRQ(ierr); 566b4319ba4SBarry Smith PetscFunctionReturn(0); 567b4319ba4SBarry Smith } 568