xref: /petsc/src/mat/utils/isltog.c (revision 2c733ebde9f595f6b5826ecb1271211dbb78ff69)
1*2c733ebdSStefano Zampini #include <petsc/private/matimpl.h>
2*2c733ebdSStefano Zampini 
3*2c733ebdSStefano Zampini PETSC_INTERN PetscErrorCode MatCreateFromISLocalToGlobalMapping(ISLocalToGlobalMapping lgmap, Mat A, PetscBool cols, PetscBool trans, MatType ptype, Mat *P)
4*2c733ebdSStefano Zampini {
5*2c733ebdSStefano Zampini   PetscBool       matfree = PETSC_FALSE;
6*2c733ebdSStefano Zampini   const PetscInt *idxs;
7*2c733ebdSStefano Zampini   PetscInt        msize, *pidxs, c = 0;
8*2c733ebdSStefano Zampini 
9*2c733ebdSStefano Zampini   PetscFunctionBegin;
10*2c733ebdSStefano Zampini   PetscValidHeaderSpecific(lgmap, IS_LTOGM_CLASSID, 1);
11*2c733ebdSStefano Zampini   PetscValidHeaderSpecific(A, MAT_CLASSID, 2);
12*2c733ebdSStefano Zampini   PetscValidType(A, 2);
13*2c733ebdSStefano Zampini   PetscValidLogicalCollectiveBool(A, cols, 3);
14*2c733ebdSStefano Zampini   PetscValidLogicalCollectiveBool(A, trans, 4);
15*2c733ebdSStefano Zampini   PetscAssertPointer(P, 6);
16*2c733ebdSStefano Zampini 
17*2c733ebdSStefano Zampini   if (!ptype) PetscCall(MatGetType(A, &ptype));
18*2c733ebdSStefano Zampini   PetscCall(PetscStrcmpAny(ptype, &matfree, MATSHELL, MATSCATTER, ""));
19*2c733ebdSStefano Zampini   PetscCall(ISLocalToGlobalMappingGetIndices(lgmap, &idxs));
20*2c733ebdSStefano Zampini   PetscCall(ISLocalToGlobalMappingGetSize(lgmap, &msize));
21*2c733ebdSStefano Zampini   PetscCall(PetscMalloc1(msize, &pidxs));
22*2c733ebdSStefano Zampini   for (PetscInt i = 0; i < msize; i++)
23*2c733ebdSStefano Zampini     if (idxs[i] >= 0) pidxs[c++] = idxs[i];
24*2c733ebdSStefano Zampini   PetscCall(ISLocalToGlobalMappingRestoreIndices(lgmap, &idxs));
25*2c733ebdSStefano Zampini   msize = c;
26*2c733ebdSStefano Zampini   if (matfree) {
27*2c733ebdSStefano Zampini     Vec        v, lv;
28*2c733ebdSStefano Zampini     VecType    vtype;
29*2c733ebdSStefano Zampini     IS         is;
30*2c733ebdSStefano Zampini     VecScatter sct;
31*2c733ebdSStefano Zampini     PetscBool  matshell;
32*2c733ebdSStefano Zampini 
33*2c733ebdSStefano Zampini     if (cols) PetscCall(MatCreateVecs(A, &v, NULL));
34*2c733ebdSStefano Zampini     else PetscCall(MatCreateVecs(A, NULL, &v));
35*2c733ebdSStefano Zampini     PetscCall(VecGetType(v, &vtype));
36*2c733ebdSStefano Zampini     PetscCall(VecCreate(PETSC_COMM_SELF, &lv));
37*2c733ebdSStefano Zampini     PetscCall(VecSetSizes(lv, msize, PETSC_DECIDE));
38*2c733ebdSStefano Zampini     PetscCall(VecSetType(lv, vtype));
39*2c733ebdSStefano Zampini     PetscCall(ISCreateGeneral(PetscObjectComm((PetscObject)A), msize, pidxs, PETSC_USE_POINTER, &is));
40*2c733ebdSStefano Zampini     if (trans) PetscCall(VecScatterCreate(lv, NULL, v, is, &sct));
41*2c733ebdSStefano Zampini     else PetscCall(VecScatterCreate(v, is, lv, NULL, &sct));
42*2c733ebdSStefano Zampini     PetscCall(MatCreateScatter(PetscObjectComm((PetscObject)A), sct, P));
43*2c733ebdSStefano Zampini     PetscCall(PetscStrcmp(ptype, MATSHELL, &matshell));
44*2c733ebdSStefano Zampini     if (matshell) {
45*2c733ebdSStefano Zampini       Mat tP;
46*2c733ebdSStefano Zampini       PetscCall(MatConvert(*P, ptype, MAT_INITIAL_MATRIX, &tP));
47*2c733ebdSStefano Zampini       PetscCall(MatDestroy(P));
48*2c733ebdSStefano Zampini       *P = tP;
49*2c733ebdSStefano Zampini     }
50*2c733ebdSStefano Zampini     PetscCall(ISDestroy(&is));
51*2c733ebdSStefano Zampini     PetscCall(VecScatterDestroy(&sct));
52*2c733ebdSStefano Zampini     PetscCall(VecDestroy(&lv));
53*2c733ebdSStefano Zampini     PetscCall(VecDestroy(&v));
54*2c733ebdSStefano Zampini   } else {
55*2c733ebdSStefano Zampini     PetscInt lar, lac, rst;
56*2c733ebdSStefano Zampini 
57*2c733ebdSStefano Zampini     PetscCall(MatGetLocalSize(A, &lar, &lac));
58*2c733ebdSStefano Zampini     PetscCall(MatCreate(PetscObjectComm((PetscObject)A), P));
59*2c733ebdSStefano Zampini     PetscCall(MatSetType(*P, ptype));
60*2c733ebdSStefano Zampini     PetscCall(MatSetSizes(*P, msize, cols ? lac : lar, PETSC_DECIDE, PETSC_DECIDE));
61*2c733ebdSStefano Zampini     PetscCall(MatSeqAIJSetPreallocation(*P, 1, NULL));
62*2c733ebdSStefano Zampini     PetscCall(MatMPIAIJSetPreallocation(*P, 1, NULL, 1, NULL));
63*2c733ebdSStefano Zampini #if defined(PETSC_HAVE_HYPRE)
64*2c733ebdSStefano Zampini     PetscCall(MatHYPRESetPreallocation(*P, 1, NULL, 1, NULL));
65*2c733ebdSStefano Zampini #endif
66*2c733ebdSStefano Zampini     PetscCall(MatGetOwnershipRange(*P, &rst, NULL));
67*2c733ebdSStefano Zampini     for (PetscInt i = 0; i < msize; i++) { PetscCall(MatSetValue(*P, i + rst, pidxs[i], 1.0, INSERT_VALUES)); }
68*2c733ebdSStefano Zampini     PetscCall(MatAssemblyBegin(*P, MAT_FINAL_ASSEMBLY));
69*2c733ebdSStefano Zampini     PetscCall(MatAssemblyEnd(*P, MAT_FINAL_ASSEMBLY));
70*2c733ebdSStefano Zampini     if (trans) {
71*2c733ebdSStefano Zampini       Mat tP;
72*2c733ebdSStefano Zampini       PetscCall(MatTranspose(*P, MAT_INITIAL_MATRIX, &tP));
73*2c733ebdSStefano Zampini       PetscCall(MatDestroy(P));
74*2c733ebdSStefano Zampini       *P = tP;
75*2c733ebdSStefano Zampini     }
76*2c733ebdSStefano Zampini   }
77*2c733ebdSStefano Zampini   PetscCall(PetscFree(pidxs));
78*2c733ebdSStefano Zampini   PetscFunctionReturn(PETSC_SUCCESS);
79*2c733ebdSStefano Zampini }
80