xref: /petsc/src/dm/impls/moab/dmmbmat.cxx (revision ccb4e88a40f0b86eaeca07ff64c64e4de2fae686)
1 #include <petsc/private/dmmbimpl.h> /*I  "petscdmmoab.h"   I*/
2 #include <petsc/private/vecimpl.h>
3 
4 #include <petscdmmoab.h>
5 #include <MBTagConventions.hpp>
6 #include <moab/NestedRefine.hpp>
7 
8 PETSC_EXTERN PetscErrorCode DMMoab_Compute_NNZ_From_Connectivity(DM, PetscInt*, PetscInt*, PetscInt*, PetscInt*, PetscBool);
9 
10 PETSC_EXTERN PetscErrorCode DMCreateMatrix_Moab(DM dm, Mat *J)
11 {
12   PetscErrorCode  ierr;
13   PetscInt        innz = 0, ionz = 0, nlsiz;
14   DM_Moab         *dmmoab = (DM_Moab*)dm->data;
15   PetscInt        *nnz = 0, *onz = 0;
16   char            *tmp = 0;
17   Mat             A;
18   MatType         mtype;
19 
20   PetscFunctionBegin;
21   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
22   PetscValidPointer(J, 3);
23 
24   /* next, need to allocate the non-zero arrays to enable pre-allocation */
25   mtype = dm->mattype;
26   ierr = PetscStrstr(mtype, MATBAIJ, &tmp);CHKERRQ(ierr);
27   nlsiz = (tmp ? dmmoab->nloc : dmmoab->nloc * dmmoab->numFields);
28 
29   /* allocate the nnz, onz arrays based on block size and local nodes */
30   ierr = PetscCalloc2(nlsiz, &nnz, nlsiz, &onz);CHKERRQ(ierr);
31 
32   /* compute the nonzero pattern based on MOAB connectivity data for local elements */
33   ierr = DMMoab_Compute_NNZ_From_Connectivity(dm, &innz, nnz, &ionz, onz, (tmp ? PETSC_TRUE : PETSC_FALSE));CHKERRQ(ierr);
34 
35   /* create the Matrix and set its type as specified by user */
36   ierr = MatCreate((((PetscObject)dm)->comm), &A);CHKERRQ(ierr);
37   ierr = MatSetSizes(A, dmmoab->nloc * dmmoab->numFields, dmmoab->nloc * dmmoab->numFields, PETSC_DETERMINE, PETSC_DETERMINE);CHKERRQ(ierr);
38   ierr = MatSetType(A, mtype);CHKERRQ(ierr);
39   ierr = MatSetBlockSize(A, dmmoab->bs);CHKERRQ(ierr);
40   ierr = MatSetDM(A, dm);CHKERRQ(ierr); /* set DM reference */
41   ierr = MatSetFromOptions(A);CHKERRQ(ierr);
42 
43   if (!dmmoab->ltog_map) SETERRQ((((PetscObject)dm)->comm), PETSC_ERR_ORDER, "Cannot create a DMMoab Mat without calling DMSetUp first.");
44   ierr = MatSetLocalToGlobalMapping(A, dmmoab->ltog_map, dmmoab->ltog_map);CHKERRQ(ierr);
45 
46   /* set preallocation based on different supported Mat types */
47   ierr = MatSeqAIJSetPreallocation(A, innz, nnz);CHKERRQ(ierr);
48   ierr = MatMPIAIJSetPreallocation(A, innz, nnz, ionz, onz);CHKERRQ(ierr);
49   ierr = MatSeqBAIJSetPreallocation(A, dmmoab->bs, innz, nnz);CHKERRQ(ierr);
50   ierr = MatMPIBAIJSetPreallocation(A, dmmoab->bs, innz, nnz, ionz, onz);CHKERRQ(ierr);
51 
52   /* clean up temporary memory */
53   ierr = PetscFree2(nnz, onz);CHKERRQ(ierr);
54 
55   /* set up internal matrix data-structures */
56   ierr = MatSetUp(A);CHKERRQ(ierr);
57 
58   /* MatSetOption(A, MAT_NEW_NONZERO_ALLOCATION_ERR, PETSC_FALSE); */
59 
60   *J = A;
61   PetscFunctionReturn(0);
62 }
63 
64 PETSC_EXTERN PetscErrorCode DMMoab_Compute_NNZ_From_Connectivity(DM dm, PetscInt* innz, PetscInt* nnz, PetscInt* ionz, PetscInt* onz, PetscBool isbaij)
65 {
66   PetscInt        i, f, nloc, vpere, bs, n_nnz, n_onz, ivtx = 0;
67   PetscInt        ibs, jbs, inbsize, iobsize, nfields, nlsiz;
68   DM_Moab         *dmmoab = (DM_Moab*)dm->data;
69   moab::Range     found;
70   std::vector<moab::EntityHandle> adjs, storage;
71   PetscBool isinterlaced;
72   moab::EntityHandle vtx;
73   moab::ErrorCode merr;
74 
75   PetscFunctionBegin;
76   bs = dmmoab->bs;
77   nloc = dmmoab->nloc;
78   nfields = dmmoab->numFields;
79   isinterlaced = (isbaij || bs == nfields ? PETSC_TRUE : PETSC_FALSE);
80   nlsiz = (isinterlaced ? nloc : nloc * nfields);
81 
82   /* loop over the locally owned vertices and figure out the NNZ pattern using connectivity information */
83   for (moab::Range::const_iterator iter = dmmoab->vowned->begin(); iter != dmmoab->vowned->end(); iter++, ivtx++) {
84 
85     vtx = *iter;
86     /* Get adjacency information for current vertex - i.e., all elements of dimension (dim) that connects
87        to the current vertex. We can then decipher if a vertex is ghosted or not and compute the
88        non-zero pattern accordingly. */
89     adjs.clear();
90     if (dmmoab->hlevel && (dmmoab->pcomm->size() == 1)) {
91       merr = dmmoab->hierarchy->get_adjacencies(vtx, dmmoab->dim, adjs); MBERRNM(merr);
92     }
93     else {
94       merr = dmmoab->mbiface->get_adjacencies(&vtx, 1, dmmoab->dim, true, adjs, moab::Interface::UNION); MBERRNM(merr);
95     }
96 
97     /* reset counters */
98     n_nnz = n_onz = 0;
99     found.clear();
100 
101     /* loop over vertices and update the number of connectivity */
102     for (unsigned jter = 0; jter < adjs.size(); ++jter) {
103 
104       /* Get connectivity information in canonical ordering for the local element */
105       const moab::EntityHandle *connect;
106       std::vector<moab::EntityHandle> cconnect;
107       merr = dmmoab->mbiface->get_connectivity(adjs[jter], connect, vpere, false, &storage); MBERRNM(merr);
108 
109       /* loop over each element connected to the adjacent vertex and update as needed */
110       for (i = 0; i < vpere; ++i) {
111         /* find the truly user-expected layer of ghosted entities to decipher NNZ pattern */
112         if (connect[i] == vtx || found.find(connect[i]) != found.end()) continue; /* make sure we don't double count shared vertices */
113         if (dmmoab->vghost->find(connect[i]) != dmmoab->vghost->end()) n_onz++; /* update out-of-proc onz */
114         else n_nnz++; /* else local vertex */
115         found.insert(connect[i]);
116       }
117     }
118     storage.clear();
119 
120     if (isbaij) {
121       nnz[ivtx] = n_nnz;  /* leave out self to avoid repeats -> node shared by multiple elements */
122       if (onz) {
123         onz[ivtx] = n_onz; /* add ghost non-owned nodes */
124       }
125     }
126     else { /* AIJ matrices */
127       if (!isinterlaced) {
128         for (f = 0; f < nfields; f++) {
129           nnz[f * nloc + ivtx] = n_nnz; /* leave out self to avoid repeats -> node shared by multiple elements */
130           if (onz)
131             onz[f * nloc + ivtx] = n_onz; /* add ghost non-owned nodes */
132         }
133       }
134       else {
135         for (f = 0; f < nfields; f++) {
136           nnz[nfields * ivtx + f] = n_nnz; /* leave out self to avoid repeats -> node shared by multiple elements */
137           if (onz)
138             onz[nfields * ivtx + f] = n_onz; /* add ghost non-owned nodes */
139         }
140       }
141     }
142   }
143 
144   for (i = 0; i < nlsiz; i++)
145     nnz[i] += 1; /* self count the node */
146 
147   for (ivtx = 0; ivtx < nloc; ivtx++) {
148     if (!isbaij) {
149       for (ibs = 0; ibs < nfields; ibs++) {
150         if (dmmoab->dfill) {  /* first address the diagonal block */
151           /* just add up the ints -- easier/faster rather than branching based on "1" */
152           for (jbs = 0, inbsize = 0; jbs < nfields; jbs++)
153             inbsize += dmmoab->dfill[ibs * nfields + jbs];
154         }
155         else inbsize = nfields; /* dense coupling since user didn't specify the component fill explicitly */
156         if (isinterlaced) nnz[ivtx * nfields + ibs] *= inbsize;
157         else nnz[ibs * nloc + ivtx] *= inbsize;
158 
159         if (onz) {
160           if (dmmoab->ofill) {  /* next address the off-diagonal block */
161             /* just add up the ints -- easier/faster rather than branching based on "1" */
162             for (jbs = 0, iobsize = 0; jbs < nfields; jbs++)
163               iobsize += dmmoab->dfill[ibs * nfields + jbs];
164           }
165           else iobsize = nfields; /* dense coupling since user didn't specify the component fill explicitly */
166           if (isinterlaced) onz[ivtx * nfields + ibs] *= iobsize;
167           else onz[ibs * nloc + ivtx] *= iobsize;
168         }
169       }
170     }
171     else {
172       /* check if we got overzealous in our nnz and onz computations */
173       nnz[ivtx] = (nnz[ivtx] > dmmoab->nloc ? dmmoab->nloc : nnz[ivtx]);
174       if (onz) onz[ivtx] = (onz[ivtx] > dmmoab->nloc ? dmmoab->nloc : onz[ivtx]);
175     }
176   }
177   /* update innz and ionz based on local maxima */
178   if (innz || ionz) {
179     if (innz) *innz = 0;
180     if (ionz) *ionz = 0;
181     for (i = 0; i < nlsiz; i++) {
182       if (innz && (nnz[i] > *innz)) *innz = nnz[i];
183       if ((ionz && onz) && (onz[i] > *ionz)) *ionz = onz[i];
184     }
185   }
186   PetscFunctionReturn(0);
187 }
188 
189 static PetscErrorCode DMMoabSetBlockFills_Private(PetscInt w, const PetscInt *fill, PetscInt **rfill)
190 {
191   PetscErrorCode ierr;
192   PetscInt       i, j, *ifill;
193 
194   PetscFunctionBegin;
195   if (!fill) PetscFunctionReturn(0);
196   ierr  = PetscMalloc1(w * w, &ifill);CHKERRQ(ierr);
197   for (i = 0; i < w; i++) {
198     for (j = 0; j < w; j++)
199       ifill[i * w + j] = fill[i * w + j];
200   }
201 
202   *rfill = ifill;
203   PetscFunctionReturn(0);
204 }
205 
206 /*@C
207     DMMoabSetBlockFills - Sets the fill pattern in each block for a multi-component problem
208     of the matrix returned by DMCreateMatrix().
209 
210     Logically Collective on da
211 
212     Input Parameters:
213 +   dm - the DMMoab object
214 .   dfill - the fill pattern in the diagonal block (may be NULL, means use dense block)
215 -   ofill - the fill pattern in the off-diagonal blocks
216 
217     Level: developer
218 
219     Notes:
220     This only makes sense when you are doing multicomponent problems but using the
221        MPIAIJ matrix format
222 
223            The format for dfill and ofill is a 2 dimensional dof by dof matrix with 1 entries
224        representing coupling and 0 entries for missing coupling. For example
225 $             dfill[9] = {1, 0, 0,
226 $                         1, 1, 0,
227 $                         0, 1, 1}
228        means that row 0 is coupled with only itself in the diagonal block, row 1 is coupled with
229        itself and row 0 (in the diagonal block) and row 2 is coupled with itself and row 1 (in the
230        diagonal block).
231 
232      DMDASetGetMatrix() allows you to provide general code for those more complicated nonzero patterns then
233      can be represented in the dfill, ofill format
234 
235    Contributed by Glenn Hammond
236 
237 .seealso DMCreateMatrix(), DMDASetGetMatrix(), DMSetMatrixPreallocateOnly()
238 
239 @*/
240 PetscErrorCode  DMMoabSetBlockFills(DM dm, const PetscInt *dfill, const PetscInt *ofill)
241 {
242   DM_Moab       *dmmoab = (DM_Moab*)dm->data;
243   PetscErrorCode ierr;
244 
245   PetscFunctionBegin;
246   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
247   ierr = DMMoabSetBlockFills_Private(dmmoab->numFields, dfill, &dmmoab->dfill);CHKERRQ(ierr);
248   ierr = DMMoabSetBlockFills_Private(dmmoab->numFields, ofill, &dmmoab->ofill);CHKERRQ(ierr);
249   PetscFunctionReturn(0);
250 }
251