xref: /petsc/src/ksp/pc/impls/gamg/util.c (revision feefa0e191a340680bb02e1467a36facdcb0b150)
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:
13*feefa0e1SJacob Faibussowitsch + Gmat    - MPIAIJ matrix for scatters
146618991cSMark Adams . data_sz - number of data terms per node (# cols in output)
15*feefa0e1SJacob Faibussowitsch    - data_in[nloc*data_sz] - column oriented local data
16bae903cbSmarkadams4 
176618991cSMark Adams   Output Parameter:
18*feefa0e1SJacob Faibussowitsch + a_stride - number of rows of output (locals+ghosts)
19*feefa0e1SJacob Faibussowitsch    - a_data_out[stride*data_sz] - output data with ghosts
20bae903cbSmarkadams4 
216618991cSMark Adams */
22d71ae5a4SJacob Faibussowitsch PetscErrorCode PCGAMGGetDataWithGhosts(Mat Gmat, PetscInt data_sz, PetscReal data_in[], PetscInt *a_stride, PetscReal **a_data_out)
23d71ae5a4SJacob Faibussowitsch {
246618991cSMark Adams   Vec          tmp_crds;
256618991cSMark Adams   Mat_MPIAIJ  *mpimat = (Mat_MPIAIJ *)Gmat->data;
266618991cSMark Adams   PetscInt     nnodes, num_ghosts, dir, kk, jj, my0, Iend, nloc;
276618991cSMark Adams   PetscScalar *data_arr;
286618991cSMark Adams   PetscReal   *datas;
296618991cSMark Adams   PetscBool    isMPIAIJ;
306618991cSMark Adams 
316618991cSMark Adams   PetscFunctionBegin;
329566063dSJacob Faibussowitsch   PetscCall(PetscObjectBaseTypeCompare((PetscObject)Gmat, MATMPIAIJ, &isMPIAIJ));
339566063dSJacob Faibussowitsch   PetscCall(MatGetOwnershipRange(Gmat, &my0, &Iend));
346618991cSMark Adams   nloc = Iend - my0;
359566063dSJacob Faibussowitsch   PetscCall(VecGetLocalSize(mpimat->lvec, &num_ghosts));
366618991cSMark Adams   nnodes    = num_ghosts + nloc;
376618991cSMark Adams   *a_stride = nnodes;
389566063dSJacob Faibussowitsch   PetscCall(MatCreateVecs(Gmat, &tmp_crds, NULL));
396618991cSMark Adams 
409566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(data_sz * nnodes, &datas));
416618991cSMark Adams   for (dir = 0; dir < data_sz; dir++) {
426618991cSMark Adams     /* set local, and global */
436618991cSMark Adams     for (kk = 0; kk < nloc; kk++) {
446618991cSMark Adams       PetscInt    gid          = my0 + kk;
456618991cSMark Adams       PetscScalar crd          = (PetscScalar)data_in[dir * nloc + kk]; /* col oriented */
46bae903cbSmarkadams4       datas[dir * nnodes + kk] = PetscRealPart(crd);                    // get local part now
476618991cSMark Adams 
489566063dSJacob Faibussowitsch       PetscCall(VecSetValues(tmp_crds, 1, &gid, &crd, INSERT_VALUES));
496618991cSMark Adams     }
509566063dSJacob Faibussowitsch     PetscCall(VecAssemblyBegin(tmp_crds));
519566063dSJacob Faibussowitsch     PetscCall(VecAssemblyEnd(tmp_crds));
52bae903cbSmarkadams4     /* scatter / gather ghost data and add to end of output data */
539566063dSJacob Faibussowitsch     PetscCall(VecScatterBegin(mpimat->Mvctx, tmp_crds, mpimat->lvec, INSERT_VALUES, SCATTER_FORWARD));
549566063dSJacob Faibussowitsch     PetscCall(VecScatterEnd(mpimat->Mvctx, tmp_crds, mpimat->lvec, INSERT_VALUES, SCATTER_FORWARD));
559566063dSJacob Faibussowitsch     PetscCall(VecGetArray(mpimat->lvec, &data_arr));
566618991cSMark Adams     for (kk = nloc, jj = 0; jj < num_ghosts; kk++, jj++) datas[dir * nnodes + kk] = PetscRealPart(data_arr[jj]);
579566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray(mpimat->lvec, &data_arr));
586618991cSMark Adams   }
599566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&tmp_crds));
606618991cSMark Adams   *a_data_out = datas;
613ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
626618991cSMark Adams }
636618991cSMark Adams 
64d71ae5a4SJacob Faibussowitsch PetscErrorCode PCGAMGHashTableCreate(PetscInt a_size, PCGAMGHashTable *a_tab)
65d71ae5a4SJacob Faibussowitsch {
666618991cSMark Adams   PetscInt kk;
676618991cSMark Adams 
686618991cSMark Adams   PetscFunctionBegin;
696618991cSMark Adams   a_tab->size = a_size;
709566063dSJacob Faibussowitsch   PetscCall(PetscMalloc2(a_size, &a_tab->table, a_size, &a_tab->data));
716618991cSMark Adams   for (kk = 0; kk < a_size; kk++) a_tab->table[kk] = -1;
723ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
736618991cSMark Adams }
746618991cSMark Adams 
75d71ae5a4SJacob Faibussowitsch PetscErrorCode PCGAMGHashTableDestroy(PCGAMGHashTable *a_tab)
76d71ae5a4SJacob Faibussowitsch {
776618991cSMark Adams   PetscFunctionBegin;
789566063dSJacob Faibussowitsch   PetscCall(PetscFree2(a_tab->table, a_tab->data));
793ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
806618991cSMark Adams }
816618991cSMark Adams 
82d71ae5a4SJacob Faibussowitsch PetscErrorCode PCGAMGHashTableAdd(PCGAMGHashTable *a_tab, PetscInt a_key, PetscInt a_data)
83d71ae5a4SJacob Faibussowitsch {
846618991cSMark Adams   PetscInt kk, idx;
856618991cSMark Adams 
866618991cSMark Adams   PetscFunctionBegin;
8763a3b9bcSJacob Faibussowitsch   PetscCheck(a_key >= 0, PETSC_COMM_SELF, PETSC_ERR_USER, "Negative key %" PetscInt_FMT ".", a_key);
888f3cd775SBarry Smith   for (kk = 0, idx = GAMG_HASH(a_key); kk < a_tab->size; kk++, idx = (idx == (a_tab->size - 1)) ? 0 : idx + 1) {
896618991cSMark Adams     if (a_tab->table[idx] == a_key) {
906618991cSMark Adams       /* exists */
916618991cSMark Adams       a_tab->data[idx] = a_data;
926618991cSMark Adams       break;
936618991cSMark Adams     } else if (a_tab->table[idx] == -1) {
946618991cSMark Adams       /* add */
956618991cSMark Adams       a_tab->table[idx] = a_key;
966618991cSMark Adams       a_tab->data[idx]  = a_data;
976618991cSMark Adams       break;
986618991cSMark Adams     }
996618991cSMark Adams   }
1006618991cSMark Adams   if (kk == a_tab->size) {
1016618991cSMark Adams     /* this is not to efficient, waiting until completely full */
1026618991cSMark Adams     PetscInt oldsize = a_tab->size, new_size = 2 * a_tab->size + 5, *oldtable = a_tab->table, *olddata = a_tab->data;
1036618991cSMark Adams 
1046618991cSMark Adams     a_tab->size = new_size;
1059566063dSJacob Faibussowitsch     PetscCall(PetscMalloc2(a_tab->size, &a_tab->table, a_tab->size, &a_tab->data));
1066618991cSMark Adams     for (kk = 0; kk < a_tab->size; kk++) a_tab->table[kk] = -1;
1076618991cSMark Adams     for (kk = 0; kk < oldsize; kk++) {
10848a46eb9SPierre Jolivet       if (oldtable[kk] != -1) PetscCall(PCGAMGHashTableAdd(a_tab, oldtable[kk], olddata[kk]));
1096618991cSMark Adams     }
1109566063dSJacob Faibussowitsch     PetscCall(PetscFree2(oldtable, olddata));
1119566063dSJacob Faibussowitsch     PetscCall(PCGAMGHashTableAdd(a_tab, a_key, a_data));
1126618991cSMark Adams   }
1133ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1146618991cSMark Adams }
115