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 804c3f3b8SBarry Smith // PetscClangLinter pragma disable: -fdoc-sowing-chars 96618991cSMark Adams /* 10bae903cbSmarkadams4 PCGAMGGetDataWithGhosts - Get array of local + ghost data with local data 11bae903cbSmarkadams4 hacks into Mat MPIAIJ so this must have size > 1 126618991cSMark Adams 1304c3f3b8SBarry Smith Input Parameters: 14feefa0e1SJacob Faibussowitsch + Gmat - MPIAIJ matrix for scatters 156618991cSMark Adams . data_sz - number of data terms per node (# cols in output) 1604c3f3b8SBarry Smith - data_in - column oriented local data of size nloc*data_sz 17bae903cbSmarkadams4 1804c3f3b8SBarry Smith Output Parameters: 19feefa0e1SJacob Faibussowitsch + a_stride - number of rows of output (locals+ghosts) 2004c3f3b8SBarry Smith - a_data_out - output data with ghosts of size stride*data_sz 21bae903cbSmarkadams4 226618991cSMark Adams */ 23d71ae5a4SJacob Faibussowitsch PetscErrorCode PCGAMGGetDataWithGhosts(Mat Gmat, PetscInt data_sz, PetscReal data_in[], PetscInt *a_stride, PetscReal **a_data_out) 24d71ae5a4SJacob Faibussowitsch { 256618991cSMark Adams Vec tmp_crds; 26*9e9caa57Smarkadams4 Mat_MPIAIJ *mpimat; 276618991cSMark Adams PetscInt nnodes, num_ghosts, dir, kk, jj, my0, Iend, nloc; 286618991cSMark Adams PetscScalar *data_arr; 296618991cSMark Adams PetscReal *datas; 306618991cSMark Adams PetscBool isMPIAIJ; 316618991cSMark Adams 326618991cSMark Adams PetscFunctionBegin; 33*9e9caa57Smarkadams4 PetscValidHeaderSpecific(Gmat, MAT_CLASSID, 1); 34*9e9caa57Smarkadams4 mpimat = (Mat_MPIAIJ *)Gmat->data; 359566063dSJacob Faibussowitsch PetscCall(PetscObjectBaseTypeCompare((PetscObject)Gmat, MATMPIAIJ, &isMPIAIJ)); 369566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(Gmat, &my0, &Iend)); 376618991cSMark Adams nloc = Iend - my0; 389566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(mpimat->lvec, &num_ghosts)); 396618991cSMark Adams nnodes = num_ghosts + nloc; 406618991cSMark Adams *a_stride = nnodes; 419566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(Gmat, &tmp_crds, NULL)); 426618991cSMark Adams 439566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(data_sz * nnodes, &datas)); 446618991cSMark Adams for (dir = 0; dir < data_sz; dir++) { 456618991cSMark Adams /* set local, and global */ 466618991cSMark Adams for (kk = 0; kk < nloc; kk++) { 476618991cSMark Adams PetscInt gid = my0 + kk; 486618991cSMark Adams PetscScalar crd = (PetscScalar)data_in[dir * nloc + kk]; /* col oriented */ 49bae903cbSmarkadams4 datas[dir * nnodes + kk] = PetscRealPart(crd); // get local part now 506618991cSMark Adams 519566063dSJacob Faibussowitsch PetscCall(VecSetValues(tmp_crds, 1, &gid, &crd, INSERT_VALUES)); 526618991cSMark Adams } 539566063dSJacob Faibussowitsch PetscCall(VecAssemblyBegin(tmp_crds)); 549566063dSJacob Faibussowitsch PetscCall(VecAssemblyEnd(tmp_crds)); 55bae903cbSmarkadams4 /* scatter / gather ghost data and add to end of output data */ 569566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(mpimat->Mvctx, tmp_crds, mpimat->lvec, INSERT_VALUES, SCATTER_FORWARD)); 579566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(mpimat->Mvctx, tmp_crds, mpimat->lvec, INSERT_VALUES, SCATTER_FORWARD)); 589566063dSJacob Faibussowitsch PetscCall(VecGetArray(mpimat->lvec, &data_arr)); 596618991cSMark Adams for (kk = nloc, jj = 0; jj < num_ghosts; kk++, jj++) datas[dir * nnodes + kk] = PetscRealPart(data_arr[jj]); 609566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(mpimat->lvec, &data_arr)); 616618991cSMark Adams } 629566063dSJacob Faibussowitsch PetscCall(VecDestroy(&tmp_crds)); 636618991cSMark Adams *a_data_out = datas; 643ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 656618991cSMark Adams } 666618991cSMark Adams 67d71ae5a4SJacob Faibussowitsch PetscErrorCode PCGAMGHashTableCreate(PetscInt a_size, PCGAMGHashTable *a_tab) 68d71ae5a4SJacob Faibussowitsch { 696618991cSMark Adams PetscInt kk; 706618991cSMark Adams 716618991cSMark Adams PetscFunctionBegin; 726618991cSMark Adams a_tab->size = a_size; 739566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(a_size, &a_tab->table, a_size, &a_tab->data)); 746618991cSMark Adams for (kk = 0; kk < a_size; kk++) a_tab->table[kk] = -1; 753ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 766618991cSMark Adams } 776618991cSMark Adams 78d71ae5a4SJacob Faibussowitsch PetscErrorCode PCGAMGHashTableDestroy(PCGAMGHashTable *a_tab) 79d71ae5a4SJacob Faibussowitsch { 806618991cSMark Adams PetscFunctionBegin; 819566063dSJacob Faibussowitsch PetscCall(PetscFree2(a_tab->table, a_tab->data)); 823ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 836618991cSMark Adams } 846618991cSMark Adams 85d71ae5a4SJacob Faibussowitsch PetscErrorCode PCGAMGHashTableAdd(PCGAMGHashTable *a_tab, PetscInt a_key, PetscInt a_data) 86d71ae5a4SJacob Faibussowitsch { 876618991cSMark Adams PetscInt kk, idx; 886618991cSMark Adams 896618991cSMark Adams PetscFunctionBegin; 9063a3b9bcSJacob Faibussowitsch PetscCheck(a_key >= 0, PETSC_COMM_SELF, PETSC_ERR_USER, "Negative key %" PetscInt_FMT ".", a_key); 918f3cd775SBarry Smith for (kk = 0, idx = GAMG_HASH(a_key); kk < a_tab->size; kk++, idx = (idx == (a_tab->size - 1)) ? 0 : idx + 1) { 926618991cSMark Adams if (a_tab->table[idx] == a_key) { 936618991cSMark Adams /* exists */ 946618991cSMark Adams a_tab->data[idx] = a_data; 956618991cSMark Adams break; 966618991cSMark Adams } else if (a_tab->table[idx] == -1) { 976618991cSMark Adams /* add */ 986618991cSMark Adams a_tab->table[idx] = a_key; 996618991cSMark Adams a_tab->data[idx] = a_data; 1006618991cSMark Adams break; 1016618991cSMark Adams } 1026618991cSMark Adams } 1036618991cSMark Adams if (kk == a_tab->size) { 1046618991cSMark Adams /* this is not to efficient, waiting until completely full */ 1056618991cSMark Adams PetscInt oldsize = a_tab->size, new_size = 2 * a_tab->size + 5, *oldtable = a_tab->table, *olddata = a_tab->data; 1066618991cSMark Adams 1076618991cSMark Adams a_tab->size = new_size; 1089566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(a_tab->size, &a_tab->table, a_tab->size, &a_tab->data)); 1096618991cSMark Adams for (kk = 0; kk < a_tab->size; kk++) a_tab->table[kk] = -1; 1106618991cSMark Adams for (kk = 0; kk < oldsize; kk++) { 11148a46eb9SPierre Jolivet if (oldtable[kk] != -1) PetscCall(PCGAMGHashTableAdd(a_tab, oldtable[kk], olddata[kk])); 1126618991cSMark Adams } 1139566063dSJacob Faibussowitsch PetscCall(PetscFree2(oldtable, olddata)); 1149566063dSJacob Faibussowitsch PetscCall(PCGAMGHashTableAdd(a_tab, a_key, a_data)); 1156618991cSMark Adams } 1163ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1176618991cSMark Adams } 118