16618991cSMark Adams /* 26618991cSMark Adams GAMG geometric-algebric multigrid PC - Mark Adams 2011 36618991cSMark Adams */ 46618991cSMark Adams #include <petsc/private/matimpl.h> 56618991cSMark Adams #include <../src/ksp/pc/impls/gamg/gamg.h> /*I "petscpc.h" I*/ 685cd6069SMark Adams #include <petsc/private/kspimpl.h> 76618991cSMark Adams 86618991cSMark Adams /* 9bae903cbSmarkadams4 PCGAMGGetDataWithGhosts - Get array of local + ghost data with local data 10bae903cbSmarkadams4 hacks into Mat MPIAIJ so this must have size > 1 116618991cSMark Adams 126618991cSMark Adams Input Parameter: 13bae903cbSmarkadams4 . Gmat - MPIAIJ matrix for scatters 146618991cSMark Adams . data_sz - number of data terms per node (# cols in output) 15bae903cbSmarkadams4 . data_in[nloc*data_sz] - column oriented local data 16bae903cbSmarkadams4 176618991cSMark Adams Output Parameter: 18bae903cbSmarkadams4 . a_stride - number of rows of output (locals+ghosts) 196618991cSMark Adams . a_data_out[stride*data_sz] - output data with ghosts 20bae903cbSmarkadams4 216618991cSMark Adams */ 229371c9d4SSatish Balay PetscErrorCode PCGAMGGetDataWithGhosts(Mat Gmat, PetscInt data_sz, PetscReal data_in[], PetscInt *a_stride, PetscReal **a_data_out) { 236618991cSMark Adams Vec tmp_crds; 246618991cSMark Adams Mat_MPIAIJ *mpimat = (Mat_MPIAIJ *)Gmat->data; 256618991cSMark Adams PetscInt nnodes, num_ghosts, dir, kk, jj, my0, Iend, nloc; 266618991cSMark Adams PetscScalar *data_arr; 276618991cSMark Adams PetscReal *datas; 286618991cSMark Adams PetscBool isMPIAIJ; 296618991cSMark Adams 306618991cSMark Adams PetscFunctionBegin; 319566063dSJacob Faibussowitsch PetscCall(PetscObjectBaseTypeCompare((PetscObject)Gmat, MATMPIAIJ, &isMPIAIJ)); 329566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(Gmat, &my0, &Iend)); 336618991cSMark Adams nloc = Iend - my0; 349566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(mpimat->lvec, &num_ghosts)); 356618991cSMark Adams nnodes = num_ghosts + nloc; 366618991cSMark Adams *a_stride = nnodes; 379566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(Gmat, &tmp_crds, NULL)); 386618991cSMark Adams 399566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(data_sz * nnodes, &datas)); 406618991cSMark Adams for (dir = 0; dir < data_sz; dir++) { 416618991cSMark Adams /* set local, and global */ 426618991cSMark Adams for (kk = 0; kk < nloc; kk++) { 436618991cSMark Adams PetscInt gid = my0 + kk; 446618991cSMark Adams PetscScalar crd = (PetscScalar)data_in[dir * nloc + kk]; /* col oriented */ 45bae903cbSmarkadams4 datas[dir * nnodes + kk] = PetscRealPart(crd); // get local part now 466618991cSMark Adams 479566063dSJacob Faibussowitsch PetscCall(VecSetValues(tmp_crds, 1, &gid, &crd, INSERT_VALUES)); 486618991cSMark Adams } 499566063dSJacob Faibussowitsch PetscCall(VecAssemblyBegin(tmp_crds)); 509566063dSJacob Faibussowitsch PetscCall(VecAssemblyEnd(tmp_crds)); 51bae903cbSmarkadams4 /* scatter / gather ghost data and add to end of output data */ 529566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(mpimat->Mvctx, tmp_crds, mpimat->lvec, INSERT_VALUES, SCATTER_FORWARD)); 539566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(mpimat->Mvctx, tmp_crds, mpimat->lvec, INSERT_VALUES, SCATTER_FORWARD)); 549566063dSJacob Faibussowitsch PetscCall(VecGetArray(mpimat->lvec, &data_arr)); 556618991cSMark Adams for (kk = nloc, jj = 0; jj < num_ghosts; kk++, jj++) datas[dir * nnodes + kk] = PetscRealPart(data_arr[jj]); 569566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(mpimat->lvec, &data_arr)); 576618991cSMark Adams } 589566063dSJacob Faibussowitsch PetscCall(VecDestroy(&tmp_crds)); 596618991cSMark Adams *a_data_out = datas; 606618991cSMark Adams PetscFunctionReturn(0); 616618991cSMark Adams } 626618991cSMark Adams 639371c9d4SSatish Balay PetscErrorCode PCGAMGHashTableCreate(PetscInt a_size, PCGAMGHashTable *a_tab) { 646618991cSMark Adams PetscInt kk; 656618991cSMark Adams 666618991cSMark Adams PetscFunctionBegin; 676618991cSMark Adams a_tab->size = a_size; 689566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(a_size, &a_tab->table, a_size, &a_tab->data)); 696618991cSMark Adams for (kk = 0; kk < a_size; kk++) a_tab->table[kk] = -1; 706618991cSMark Adams PetscFunctionReturn(0); 716618991cSMark Adams } 726618991cSMark Adams 739371c9d4SSatish Balay PetscErrorCode PCGAMGHashTableDestroy(PCGAMGHashTable *a_tab) { 746618991cSMark Adams PetscFunctionBegin; 759566063dSJacob Faibussowitsch PetscCall(PetscFree2(a_tab->table, a_tab->data)); 766618991cSMark Adams PetscFunctionReturn(0); 776618991cSMark Adams } 786618991cSMark Adams 799371c9d4SSatish Balay PetscErrorCode PCGAMGHashTableAdd(PCGAMGHashTable *a_tab, PetscInt a_key, PetscInt a_data) { 806618991cSMark Adams PetscInt kk, idx; 816618991cSMark Adams 826618991cSMark Adams PetscFunctionBegin; 8363a3b9bcSJacob Faibussowitsch PetscCheck(a_key >= 0, PETSC_COMM_SELF, PETSC_ERR_USER, "Negative key %" PetscInt_FMT ".", a_key); 848f3cd775SBarry Smith for (kk = 0, idx = GAMG_HASH(a_key); kk < a_tab->size; kk++, idx = (idx == (a_tab->size - 1)) ? 0 : idx + 1) { 856618991cSMark Adams if (a_tab->table[idx] == a_key) { 866618991cSMark Adams /* exists */ 876618991cSMark Adams a_tab->data[idx] = a_data; 886618991cSMark Adams break; 896618991cSMark Adams } else if (a_tab->table[idx] == -1) { 906618991cSMark Adams /* add */ 916618991cSMark Adams a_tab->table[idx] = a_key; 926618991cSMark Adams a_tab->data[idx] = a_data; 936618991cSMark Adams break; 946618991cSMark Adams } 956618991cSMark Adams } 966618991cSMark Adams if (kk == a_tab->size) { 976618991cSMark Adams /* this is not to efficient, waiting until completely full */ 986618991cSMark Adams PetscInt oldsize = a_tab->size, new_size = 2 * a_tab->size + 5, *oldtable = a_tab->table, *olddata = a_tab->data; 996618991cSMark Adams 1006618991cSMark Adams a_tab->size = new_size; 1019566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(a_tab->size, &a_tab->table, a_tab->size, &a_tab->data)); 1026618991cSMark Adams for (kk = 0; kk < a_tab->size; kk++) a_tab->table[kk] = -1; 1036618991cSMark Adams for (kk = 0; kk < oldsize; kk++) { 104*48a46eb9SPierre Jolivet if (oldtable[kk] != -1) PetscCall(PCGAMGHashTableAdd(a_tab, oldtable[kk], olddata[kk])); 1056618991cSMark Adams } 1069566063dSJacob Faibussowitsch PetscCall(PetscFree2(oldtable, olddata)); 1079566063dSJacob Faibussowitsch PetscCall(PCGAMGHashTableAdd(a_tab, a_key, a_data)); 1086618991cSMark Adams } 1096618991cSMark Adams PetscFunctionReturn(0); 1106618991cSMark Adams } 111