xref: /petsc/src/dm/interface/dmi.c (revision e611a964e9853b74d61a56642fe9d06a6e51780f)
1 #include <petsc/private/dmimpl.h>     /*I      "petscdm.h"     I*/
2 #include <petscds.h>
3 
4 #undef __FUNCT__
5 #define __FUNCT__ "DMCreateGlobalVector_Section_Private"
6 PetscErrorCode DMCreateGlobalVector_Section_Private(DM dm,Vec *vec)
7 {
8   PetscSection   gSection;
9   PetscInt       localSize, bs, blockSize = -1, pStart, pEnd, p;
10   PetscErrorCode ierr;
11 
12   PetscFunctionBegin;
13   ierr = DMGetDefaultGlobalSection(dm, &gSection);CHKERRQ(ierr);
14   ierr = PetscSectionGetChart(gSection, &pStart, &pEnd);CHKERRQ(ierr);
15   for (p = pStart; p < pEnd; ++p) {
16     PetscInt dof, cdof;
17 
18     ierr = PetscSectionGetDof(gSection, p, &dof);CHKERRQ(ierr);
19     ierr = PetscSectionGetConstraintDof(gSection, p, &cdof);CHKERRQ(ierr);
20     if ((blockSize < 0) && (dof > 0) && (dof-cdof > 0)) blockSize = dof-cdof;
21     if ((dof > 0) && (dof-cdof != blockSize)) {
22       blockSize = 1;
23       break;
24     }
25   }
26   if (blockSize < 0) blockSize = PETSC_MAX_INT;
27   ierr = MPIU_Allreduce(&blockSize, &bs, 1, MPIU_INT, MPI_MIN, PetscObjectComm((PetscObject)dm));CHKERRQ(ierr);
28   if (blockSize == PETSC_MAX_INT) blockSize = 1; /* Everyone was empty */
29   ierr = PetscSectionGetConstrainedStorageSize(gSection, &localSize);CHKERRQ(ierr);
30   if (localSize%blockSize) SETERRQ2(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Mismatch between blocksize %d and local storage size %d", blockSize, localSize);
31   ierr = VecCreate(PetscObjectComm((PetscObject)dm), vec);CHKERRQ(ierr);
32   ierr = VecSetSizes(*vec, localSize, PETSC_DETERMINE);CHKERRQ(ierr);
33   ierr = VecSetBlockSize(*vec, bs);CHKERRQ(ierr);
34   ierr = VecSetType(*vec,dm->vectype);CHKERRQ(ierr);
35   ierr = VecSetDM(*vec, dm);CHKERRQ(ierr);
36   /* ierr = VecSetLocalToGlobalMapping(*vec, dm->ltogmap);CHKERRQ(ierr); */
37   PetscFunctionReturn(0);
38 }
39 
40 #undef __FUNCT__
41 #define __FUNCT__ "DMCreateLocalVector_Section_Private"
42 PetscErrorCode DMCreateLocalVector_Section_Private(DM dm,Vec *vec)
43 {
44   PetscSection   section;
45   PetscInt       localSize, blockSize = -1, pStart, pEnd, p;
46   PetscErrorCode ierr;
47 
48   PetscFunctionBegin;
49   ierr = DMGetDefaultSection(dm, &section);CHKERRQ(ierr);
50   ierr = PetscSectionGetChart(section, &pStart, &pEnd);CHKERRQ(ierr);
51   for (p = pStart; p < pEnd; ++p) {
52     PetscInt dof;
53 
54     ierr = PetscSectionGetDof(section, p, &dof);CHKERRQ(ierr);
55     if ((blockSize < 0) && (dof > 0)) blockSize = dof;
56     if ((dof > 0) && (dof != blockSize)) {
57       blockSize = 1;
58       break;
59     }
60   }
61   ierr = PetscSectionGetStorageSize(section, &localSize);CHKERRQ(ierr);
62   ierr = VecCreate(PETSC_COMM_SELF, vec);CHKERRQ(ierr);
63   ierr = VecSetSizes(*vec, localSize, localSize);CHKERRQ(ierr);
64   ierr = VecSetBlockSize(*vec, blockSize);CHKERRQ(ierr);
65   ierr = VecSetType(*vec,dm->vectype);CHKERRQ(ierr);
66   ierr = VecSetDM(*vec, dm);CHKERRQ(ierr);
67   PetscFunctionReturn(0);
68 }
69 
70 #undef __FUNCT__
71 #define __FUNCT__ "DMCreateSubDM_Section_Private"
72 /* This assumes that the DM has been cloned prior to the call */
73 PetscErrorCode DMCreateSubDM_Section_Private(DM dm, PetscInt numFields, PetscInt fields[], IS *is, DM *subdm)
74 {
75   PetscSection   section, sectionGlobal;
76   PetscInt      *subIndices;
77   PetscInt       subSize = 0, subOff = 0, nF, f, pStart, pEnd, p;
78   PetscErrorCode ierr;
79 
80   PetscFunctionBegin;
81   if (!numFields) PetscFunctionReturn(0);
82   ierr = DMGetDefaultSection(dm, &section);CHKERRQ(ierr);
83   ierr = DMGetDefaultGlobalSection(dm, &sectionGlobal);CHKERRQ(ierr);
84   if (!section) SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Must set default section for DM before splitting fields");
85   if (!sectionGlobal) SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Must set default global section for DM before splitting fields");
86   ierr = PetscSectionGetNumFields(section, &nF);CHKERRQ(ierr);
87   if (numFields > nF) SETERRQ2(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Number of requested fields %d greater than number of DM fields %d", numFields, nF);
88   if (is) {
89     PetscInt bs = -1, bsLocal, bsMax, bsMin;
90     ierr = PetscSectionGetChart(sectionGlobal, &pStart, &pEnd);CHKERRQ(ierr);
91     for (p = pStart; p < pEnd; ++p) {
92       PetscInt gdof, pSubSize  = 0;
93 
94       ierr = PetscSectionGetDof(sectionGlobal, p, &gdof);CHKERRQ(ierr);
95       if (gdof > 0) {
96         for (f = 0; f < numFields; ++f) {
97           PetscInt fdof, fcdof;
98 
99           ierr     = PetscSectionGetFieldDof(section, p, fields[f], &fdof);CHKERRQ(ierr);
100           ierr     = PetscSectionGetFieldConstraintDof(section, p, fields[f], &fcdof);CHKERRQ(ierr);
101           pSubSize += fdof-fcdof;
102         }
103         subSize += pSubSize;
104         if (pSubSize) {
105           if (bs < 0) {
106             bs = pSubSize;
107           } else if (bs != pSubSize) {
108             /* Layout does not admit a pointwise block size */
109             bs = 1;
110           }
111         }
112       }
113     }
114     /* Must have same blocksize on all procs (some might have no points) */
115     bsLocal = bs;
116     ierr = MPIU_Allreduce(&bsLocal, &bsMax, 1, MPIU_INT, MPI_MAX, PetscObjectComm((PetscObject)dm));CHKERRQ(ierr);
117     bsLocal = bs < 0 ? bsMax : bs;
118     ierr = MPIU_Allreduce(&bsLocal, &bsMin, 1, MPIU_INT, MPI_MIN, PetscObjectComm((PetscObject)dm));CHKERRQ(ierr);
119     if (bsMin != bsMax) {
120       bs = 1;
121     } else {
122       bs = bsMax;
123     }
124     ierr = PetscMalloc1(subSize, &subIndices);CHKERRQ(ierr);
125     for (p = pStart; p < pEnd; ++p) {
126       PetscInt gdof, goff;
127 
128       ierr = PetscSectionGetDof(sectionGlobal, p, &gdof);CHKERRQ(ierr);
129       if (gdof > 0) {
130         ierr = PetscSectionGetOffset(sectionGlobal, p, &goff);CHKERRQ(ierr);
131         for (f = 0; f < numFields; ++f) {
132           PetscInt fdof, fcdof, fc, f2, poff = 0;
133 
134           /* Can get rid of this loop by storing field information in the global section */
135           for (f2 = 0; f2 < fields[f]; ++f2) {
136             ierr  = PetscSectionGetFieldDof(section, p, f2, &fdof);CHKERRQ(ierr);
137             ierr  = PetscSectionGetFieldConstraintDof(section, p, f2, &fcdof);CHKERRQ(ierr);
138             poff += fdof-fcdof;
139           }
140           ierr = PetscSectionGetFieldDof(section, p, fields[f], &fdof);CHKERRQ(ierr);
141           ierr = PetscSectionGetFieldConstraintDof(section, p, fields[f], &fcdof);CHKERRQ(ierr);
142           for (fc = 0; fc < fdof-fcdof; ++fc, ++subOff) {
143             subIndices[subOff] = goff+poff+fc;
144           }
145         }
146       }
147     }
148     ierr = ISCreateGeneral(PetscObjectComm((PetscObject)dm), subSize, subIndices, PETSC_OWN_POINTER, is);CHKERRQ(ierr);
149     ierr = ISSetBlockSize(*is, bs);CHKERRQ(ierr);
150   }
151   if (subdm) {
152     PetscSection subsection;
153     PetscBool    haveNull = PETSC_FALSE;
154     PetscInt     f, nf = 0;
155 
156     ierr = PetscSectionCreateSubsection(section, numFields, fields, &subsection);CHKERRQ(ierr);
157     ierr = DMSetDefaultSection(*subdm, subsection);CHKERRQ(ierr);
158     ierr = PetscSectionDestroy(&subsection);CHKERRQ(ierr);
159     for (f = 0; f < numFields; ++f) {
160       (*subdm)->nullspaceConstructors[f] = dm->nullspaceConstructors[fields[f]];
161       if ((*subdm)->nullspaceConstructors[f]) {
162         haveNull = PETSC_TRUE;
163         nf       = f;
164       }
165     }
166     if (haveNull && is) {
167       MatNullSpace nullSpace;
168 
169       ierr = (*(*subdm)->nullspaceConstructors[nf])(*subdm, nf, &nullSpace);CHKERRQ(ierr);
170       ierr = PetscObjectCompose((PetscObject) *is, "nullspace", (PetscObject) nullSpace);CHKERRQ(ierr);
171       ierr = MatNullSpaceDestroy(&nullSpace);CHKERRQ(ierr);
172     }
173     if (dm->prob) {
174       PetscInt Nf;
175 
176       ierr = PetscDSGetNumFields(dm->prob, &Nf);CHKERRQ(ierr);
177       if (nF != Nf) SETERRQ2(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "The number of DM fields %d does not match the number of Section fields %d", Nf, nF);
178       ierr = DMSetNumFields(*subdm, numFields);CHKERRQ(ierr);
179       for (f = 0; f < numFields; ++f) {
180         PetscObject disc;
181 
182         ierr = DMGetField(dm, fields[f], &disc);CHKERRQ(ierr);
183         ierr = DMSetField(*subdm, f, disc);CHKERRQ(ierr);
184       }
185       if (numFields == 1 && is) {
186         PetscObject disc, space, pmat;
187 
188         ierr = DMGetField(*subdm, 0, &disc);CHKERRQ(ierr);
189         ierr = PetscObjectQuery(disc, "nullspace", &space);CHKERRQ(ierr);
190         if (space) {ierr = PetscObjectCompose((PetscObject) *is, "nullspace", space);CHKERRQ(ierr);}
191         ierr = PetscObjectQuery(disc, "nearnullspace", &space);CHKERRQ(ierr);
192         if (space) {ierr = PetscObjectCompose((PetscObject) *is, "nearnullspace", space);CHKERRQ(ierr);}
193         ierr = PetscObjectQuery(disc, "pmat", &pmat);CHKERRQ(ierr);
194         if (pmat) {ierr = PetscObjectCompose((PetscObject) *is, "pmat", pmat);CHKERRQ(ierr);}
195       }
196     }
197   }
198 #if 0
199   /* We need a way to filter the original SF for given fields:
200        - Keeping the original section around it too much I think
201        - We could keep the distributed section, and subset it
202    */
203   if (dm->sfNatural) {
204     PetscSF sfNatural;
205 
206     ierr = PetscSectionCreateSubsection(dm->originalSection, numFields, fields, &(*subdm)->originalSection);CHKERRQ(ierr);
207     ierr = DMPlexCreateGlobalToNaturalPetscSF(*subdm, &sfNatural);CHKERRQ(ierr);
208     ierr = DMPlexSetGlobalToNaturalPetscSF(*subdm, sfNatural);CHKERRQ(ierr);
209   }
210 #endif
211   PetscFunctionReturn(0);
212 }
213