12c733ebdSStefano Zampini #include <petsc/private/matimpl.h>
22c733ebdSStefano Zampini
MatCreateFromISLocalToGlobalMapping(ISLocalToGlobalMapping lgmap,Mat A,PetscBool cols,PetscBool trans,MatType ptype,Mat * P)3239f794eSPierre Jolivet PetscErrorCode MatCreateFromISLocalToGlobalMapping(ISLocalToGlobalMapping lgmap, Mat A, PetscBool cols, PetscBool trans, MatType ptype, Mat *P)
42c733ebdSStefano Zampini {
52c733ebdSStefano Zampini PetscBool matfree = PETSC_FALSE;
62c733ebdSStefano Zampini const PetscInt *idxs;
72c733ebdSStefano Zampini PetscInt msize, *pidxs, c = 0;
82c733ebdSStefano Zampini
92c733ebdSStefano Zampini PetscFunctionBegin;
102c733ebdSStefano Zampini PetscValidHeaderSpecific(lgmap, IS_LTOGM_CLASSID, 1);
112c733ebdSStefano Zampini PetscValidHeaderSpecific(A, MAT_CLASSID, 2);
122c733ebdSStefano Zampini PetscValidType(A, 2);
132c733ebdSStefano Zampini PetscValidLogicalCollectiveBool(A, cols, 3);
142c733ebdSStefano Zampini PetscValidLogicalCollectiveBool(A, trans, 4);
152c733ebdSStefano Zampini PetscAssertPointer(P, 6);
162c733ebdSStefano Zampini
172c733ebdSStefano Zampini if (!ptype) PetscCall(MatGetType(A, &ptype));
182c733ebdSStefano Zampini PetscCall(PetscStrcmpAny(ptype, &matfree, MATSHELL, MATSCATTER, ""));
192c733ebdSStefano Zampini PetscCall(ISLocalToGlobalMappingGetIndices(lgmap, &idxs));
202c733ebdSStefano Zampini PetscCall(ISLocalToGlobalMappingGetSize(lgmap, &msize));
212c733ebdSStefano Zampini PetscCall(PetscMalloc1(msize, &pidxs));
222c733ebdSStefano Zampini for (PetscInt i = 0; i < msize; i++)
232c733ebdSStefano Zampini if (idxs[i] >= 0) pidxs[c++] = idxs[i];
242c733ebdSStefano Zampini PetscCall(ISLocalToGlobalMappingRestoreIndices(lgmap, &idxs));
252c733ebdSStefano Zampini msize = c;
262c733ebdSStefano Zampini if (matfree) {
272c733ebdSStefano Zampini Vec v, lv;
282c733ebdSStefano Zampini VecType vtype;
292c733ebdSStefano Zampini IS is;
302c733ebdSStefano Zampini VecScatter sct;
312c733ebdSStefano Zampini PetscBool matshell;
322c733ebdSStefano Zampini
332c733ebdSStefano Zampini if (cols) PetscCall(MatCreateVecs(A, &v, NULL));
342c733ebdSStefano Zampini else PetscCall(MatCreateVecs(A, NULL, &v));
352c733ebdSStefano Zampini PetscCall(VecGetType(v, &vtype));
362c733ebdSStefano Zampini PetscCall(VecCreate(PETSC_COMM_SELF, &lv));
372c733ebdSStefano Zampini PetscCall(VecSetSizes(lv, msize, PETSC_DECIDE));
382c733ebdSStefano Zampini PetscCall(VecSetType(lv, vtype));
392c733ebdSStefano Zampini PetscCall(ISCreateGeneral(PetscObjectComm((PetscObject)A), msize, pidxs, PETSC_USE_POINTER, &is));
402c733ebdSStefano Zampini if (trans) PetscCall(VecScatterCreate(lv, NULL, v, is, &sct));
412c733ebdSStefano Zampini else PetscCall(VecScatterCreate(v, is, lv, NULL, &sct));
422c733ebdSStefano Zampini PetscCall(MatCreateScatter(PetscObjectComm((PetscObject)A), sct, P));
432c733ebdSStefano Zampini PetscCall(PetscStrcmp(ptype, MATSHELL, &matshell));
442c733ebdSStefano Zampini if (matshell) {
452c733ebdSStefano Zampini Mat tP;
462c733ebdSStefano Zampini PetscCall(MatConvert(*P, ptype, MAT_INITIAL_MATRIX, &tP));
472c733ebdSStefano Zampini PetscCall(MatDestroy(P));
482c733ebdSStefano Zampini *P = tP;
492c733ebdSStefano Zampini }
502c733ebdSStefano Zampini PetscCall(ISDestroy(&is));
512c733ebdSStefano Zampini PetscCall(VecScatterDestroy(&sct));
522c733ebdSStefano Zampini PetscCall(VecDestroy(&lv));
532c733ebdSStefano Zampini PetscCall(VecDestroy(&v));
542c733ebdSStefano Zampini } else {
552c733ebdSStefano Zampini PetscInt lar, lac, rst;
562c733ebdSStefano Zampini
572c733ebdSStefano Zampini PetscCall(MatGetLocalSize(A, &lar, &lac));
582c733ebdSStefano Zampini PetscCall(MatCreate(PetscObjectComm((PetscObject)A), P));
592c733ebdSStefano Zampini PetscCall(MatSetType(*P, ptype));
602c733ebdSStefano Zampini PetscCall(MatSetSizes(*P, msize, cols ? lac : lar, PETSC_DECIDE, PETSC_DECIDE));
612c733ebdSStefano Zampini PetscCall(MatSeqAIJSetPreallocation(*P, 1, NULL));
622c733ebdSStefano Zampini PetscCall(MatMPIAIJSetPreallocation(*P, 1, NULL, 1, NULL));
632c733ebdSStefano Zampini #if defined(PETSC_HAVE_HYPRE)
642c733ebdSStefano Zampini PetscCall(MatHYPRESetPreallocation(*P, 1, NULL, 1, NULL));
652c733ebdSStefano Zampini #endif
662c733ebdSStefano Zampini PetscCall(MatGetOwnershipRange(*P, &rst, NULL));
67*3a7d0413SPierre Jolivet for (PetscInt i = 0; i < msize; i++) PetscCall(MatSetValue(*P, i + rst, pidxs[i], 1.0, INSERT_VALUES));
682c733ebdSStefano Zampini PetscCall(MatAssemblyBegin(*P, MAT_FINAL_ASSEMBLY));
692c733ebdSStefano Zampini PetscCall(MatAssemblyEnd(*P, MAT_FINAL_ASSEMBLY));
702c733ebdSStefano Zampini if (trans) {
712c733ebdSStefano Zampini Mat tP;
722c733ebdSStefano Zampini PetscCall(MatTranspose(*P, MAT_INITIAL_MATRIX, &tP));
732c733ebdSStefano Zampini PetscCall(MatDestroy(P));
742c733ebdSStefano Zampini *P = tP;
752c733ebdSStefano Zampini }
762c733ebdSStefano Zampini }
772c733ebdSStefano Zampini PetscCall(PetscFree(pidxs));
782c733ebdSStefano Zampini PetscFunctionReturn(PETSC_SUCCESS);
792c733ebdSStefano Zampini }
80