xref: /petsc/src/dm/interface/dmi.c (revision 4d9407bc7dfa10a1d8b044e867a3008fab670d40)
111689aebSJed Brown #include <petsc-private/dmimpl.h>     /*I      "petscdm.h"     I*/
211689aebSJed Brown 
311689aebSJed Brown #undef __FUNCT__
411689aebSJed Brown #define __FUNCT__ "DMCreateGlobalVector_Section_Private"
511689aebSJed Brown PetscErrorCode DMCreateGlobalVector_Section_Private(DM dm,Vec *vec)
611689aebSJed Brown {
711689aebSJed Brown   PetscSection   gSection;
8cc30b04dSMatthew G Knepley   PetscInt       localSize, bs, blockSize = -1, pStart, pEnd, p;
9fad22124SMatthew G Knepley   PetscErrorCode ierr;
1011689aebSJed Brown 
1111689aebSJed Brown   PetscFunctionBegin;
1211689aebSJed Brown   ierr = DMGetDefaultGlobalSection(dm, &gSection);CHKERRQ(ierr);
13fad22124SMatthew G Knepley   ierr = PetscSectionGetChart(gSection, &pStart, &pEnd);CHKERRQ(ierr);
1411689aebSJed Brown   for (p = pStart; p < pEnd; ++p) {
1511689aebSJed Brown     PetscInt dof, cdof;
1611689aebSJed Brown 
17fad22124SMatthew G Knepley     ierr = PetscSectionGetDof(gSection, p, &dof);CHKERRQ(ierr);
18fad22124SMatthew G Knepley     ierr = PetscSectionGetConstraintDof(gSection, p, &cdof);CHKERRQ(ierr);
191f680262SMatthew G Knepley     if ((blockSize < 0) && (dof > 0) && (dof-cdof > 0)) blockSize = dof-cdof;
2011689aebSJed Brown     if ((dof > 0) && (dof-cdof != blockSize)) {
2111689aebSJed Brown       blockSize = 1;
2211689aebSJed Brown       break;
2311689aebSJed Brown     }
2411689aebSJed Brown   }
25cc30b04dSMatthew G Knepley   if (blockSize < 0) blockSize = 1;
2682f516ccSBarry Smith   ierr = MPI_Allreduce(&blockSize, &bs, 1, MPIU_INT, MPI_MIN, PetscObjectComm((PetscObject)dm));CHKERRQ(ierr);
27fad22124SMatthew G Knepley   ierr = PetscSectionGetConstrainedStorageSize(gSection, &localSize);CHKERRQ(ierr);
2882f516ccSBarry Smith   if (localSize%blockSize) SETERRQ2(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Mismatch between blocksize %d and local storage size %d", blockSize, localSize);
2982f516ccSBarry Smith   ierr = VecCreate(PetscObjectComm((PetscObject)dm), vec);CHKERRQ(ierr);
3011689aebSJed Brown   ierr = VecSetSizes(*vec, localSize, PETSC_DETERMINE);CHKERRQ(ierr);
31cc30b04dSMatthew G Knepley   ierr = VecSetBlockSize(*vec, bs);CHKERRQ(ierr);
3211689aebSJed Brown   /* ierr = VecSetType(*vec, dm->vectype);CHKERRQ(ierr); */
3311689aebSJed Brown   ierr = VecSetFromOptions(*vec);CHKERRQ(ierr);
3411689aebSJed Brown   ierr = VecSetDM(*vec, dm);CHKERRQ(ierr);
3511689aebSJed Brown   /* ierr = VecSetLocalToGlobalMapping(*vec, dm->ltogmap);CHKERRQ(ierr); */
3611689aebSJed Brown   /* ierr = VecSetLocalToGlobalMappingBlock(*vec, dm->ltogmapb);CHKERRQ(ierr); */
3711689aebSJed Brown   PetscFunctionReturn(0);
3811689aebSJed Brown }
3911689aebSJed Brown 
4011689aebSJed Brown #undef __FUNCT__
4111689aebSJed Brown #define __FUNCT__ "DMCreateLocalVector_Section_Private"
4211689aebSJed Brown PetscErrorCode DMCreateLocalVector_Section_Private(DM dm,Vec *vec)
4311689aebSJed Brown {
44fad22124SMatthew G Knepley   PetscSection   section;
4511689aebSJed Brown   PetscInt       localSize, blockSize = -1, pStart, pEnd, p;
46fad22124SMatthew G Knepley   PetscErrorCode ierr;
4711689aebSJed Brown 
4811689aebSJed Brown   PetscFunctionBegin;
49fad22124SMatthew G Knepley   ierr = DMGetDefaultSection(dm, &section);CHKERRQ(ierr);
50fad22124SMatthew G Knepley   ierr = PetscSectionGetChart(section, &pStart, &pEnd);CHKERRQ(ierr);
5111689aebSJed Brown   for (p = pStart; p < pEnd; ++p) {
5211689aebSJed Brown     PetscInt dof;
5311689aebSJed Brown 
54fad22124SMatthew G Knepley     ierr = PetscSectionGetDof(section, p, &dof);CHKERRQ(ierr);
5511689aebSJed Brown     if ((blockSize < 0) && (dof > 0)) blockSize = dof;
5611689aebSJed Brown     if ((dof > 0) && (dof != blockSize)) {
5711689aebSJed Brown       blockSize = 1;
5811689aebSJed Brown       break;
5911689aebSJed Brown     }
6011689aebSJed Brown   }
61fad22124SMatthew G Knepley   ierr = PetscSectionGetStorageSize(section, &localSize);CHKERRQ(ierr);
6211689aebSJed Brown   ierr = VecCreate(PETSC_COMM_SELF, vec);CHKERRQ(ierr);
6311689aebSJed Brown   ierr = VecSetSizes(*vec, localSize, localSize);CHKERRQ(ierr);
6411689aebSJed Brown   ierr = VecSetBlockSize(*vec, blockSize);CHKERRQ(ierr);
6511689aebSJed Brown   ierr = VecSetFromOptions(*vec);CHKERRQ(ierr);
6611689aebSJed Brown   ierr = VecSetDM(*vec, dm);CHKERRQ(ierr);
6711689aebSJed Brown   PetscFunctionReturn(0);
6811689aebSJed Brown }
69*4d9407bcSMatthew G. Knepley 
70*4d9407bcSMatthew G. Knepley #undef __FUNCT__
71*4d9407bcSMatthew G. Knepley #define __FUNCT__ "DMCreateSubDM_Section_Private"
72*4d9407bcSMatthew G. Knepley /* This assumes that the DM has been cloned prior to the call */
73*4d9407bcSMatthew G. Knepley PetscErrorCode DMCreateSubDM_Section_Private(DM dm, PetscInt numFields, PetscInt fields[], IS *is, DM *subdm)
74*4d9407bcSMatthew G. Knepley {
75*4d9407bcSMatthew G. Knepley   PetscSection   section, sectionGlobal;
76*4d9407bcSMatthew G. Knepley   PetscInt      *subIndices;
77*4d9407bcSMatthew G. Knepley   PetscInt       subSize = 0, subOff = 0, nF, f, pStart, pEnd, p;
78*4d9407bcSMatthew G. Knepley   PetscErrorCode ierr;
79*4d9407bcSMatthew G. Knepley 
80*4d9407bcSMatthew G. Knepley   PetscFunctionBegin;
81*4d9407bcSMatthew G. Knepley   if (!numFields) PetscFunctionReturn(0);
82*4d9407bcSMatthew G. Knepley   ierr = DMGetDefaultSection(dm, &section);CHKERRQ(ierr);
83*4d9407bcSMatthew G. Knepley   ierr = DMGetDefaultGlobalSection(dm, &sectionGlobal);CHKERRQ(ierr);
84*4d9407bcSMatthew G. Knepley   if (!section) SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Must set default section for DM before splitting fields");
85*4d9407bcSMatthew G. Knepley   if (!sectionGlobal) SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Must set default global section for DM before splitting fields");
86*4d9407bcSMatthew G. Knepley   ierr = PetscSectionGetNumFields(section, &nF);CHKERRQ(ierr);
87*4d9407bcSMatthew G. Knepley   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*4d9407bcSMatthew G. Knepley   if (is) {
89*4d9407bcSMatthew G. Knepley     ierr = PetscSectionGetChart(sectionGlobal, &pStart, &pEnd);CHKERRQ(ierr);
90*4d9407bcSMatthew G. Knepley     for (p = pStart; p < pEnd; ++p) {
91*4d9407bcSMatthew G. Knepley       PetscInt gdof;
92*4d9407bcSMatthew G. Knepley 
93*4d9407bcSMatthew G. Knepley       ierr = PetscSectionGetDof(sectionGlobal, p, &gdof);CHKERRQ(ierr);
94*4d9407bcSMatthew G. Knepley       if (gdof > 0) {
95*4d9407bcSMatthew G. Knepley         for (f = 0; f < numFields; ++f) {
96*4d9407bcSMatthew G. Knepley           PetscInt fdof, fcdof;
97*4d9407bcSMatthew G. Knepley 
98*4d9407bcSMatthew G. Knepley           ierr     = PetscSectionGetFieldDof(section, p, fields[f], &fdof);CHKERRQ(ierr);
99*4d9407bcSMatthew G. Knepley           ierr     = PetscSectionGetFieldConstraintDof(section, p, fields[f], &fcdof);CHKERRQ(ierr);
100*4d9407bcSMatthew G. Knepley           subSize += fdof-fcdof;
101*4d9407bcSMatthew G. Knepley         }
102*4d9407bcSMatthew G. Knepley       }
103*4d9407bcSMatthew G. Knepley     }
104*4d9407bcSMatthew G. Knepley     ierr = PetscMalloc(subSize * sizeof(PetscInt), &subIndices);CHKERRQ(ierr);
105*4d9407bcSMatthew G. Knepley     for (p = pStart; p < pEnd; ++p) {
106*4d9407bcSMatthew G. Knepley       PetscInt gdof, goff;
107*4d9407bcSMatthew G. Knepley 
108*4d9407bcSMatthew G. Knepley       ierr = PetscSectionGetDof(sectionGlobal, p, &gdof);CHKERRQ(ierr);
109*4d9407bcSMatthew G. Knepley       if (gdof > 0) {
110*4d9407bcSMatthew G. Knepley         ierr = PetscSectionGetOffset(sectionGlobal, p, &goff);CHKERRQ(ierr);
111*4d9407bcSMatthew G. Knepley         for (f = 0; f < numFields; ++f) {
112*4d9407bcSMatthew G. Knepley           PetscInt fdof, fcdof, fc, f2, poff = 0;
113*4d9407bcSMatthew G. Knepley 
114*4d9407bcSMatthew G. Knepley           /* Can get rid of this loop by storing field information in the global section */
115*4d9407bcSMatthew G. Knepley           for (f2 = 0; f2 < fields[f]; ++f2) {
116*4d9407bcSMatthew G. Knepley             ierr  = PetscSectionGetFieldDof(section, p, f2, &fdof);CHKERRQ(ierr);
117*4d9407bcSMatthew G. Knepley             ierr  = PetscSectionGetFieldConstraintDof(section, p, f2, &fcdof);CHKERRQ(ierr);
118*4d9407bcSMatthew G. Knepley             poff += fdof-fcdof;
119*4d9407bcSMatthew G. Knepley           }
120*4d9407bcSMatthew G. Knepley           ierr = PetscSectionGetFieldDof(section, p, fields[f], &fdof);CHKERRQ(ierr);
121*4d9407bcSMatthew G. Knepley           ierr = PetscSectionGetFieldConstraintDof(section, p, fields[f], &fcdof);CHKERRQ(ierr);
122*4d9407bcSMatthew G. Knepley           for (fc = 0; fc < fdof-fcdof; ++fc, ++subOff) {
123*4d9407bcSMatthew G. Knepley             subIndices[subOff] = goff+poff+fc;
124*4d9407bcSMatthew G. Knepley           }
125*4d9407bcSMatthew G. Knepley         }
126*4d9407bcSMatthew G. Knepley       }
127*4d9407bcSMatthew G. Knepley     }
128*4d9407bcSMatthew G. Knepley     ierr = ISCreateGeneral(PetscObjectComm((PetscObject)dm), subSize, subIndices, PETSC_OWN_POINTER, is);CHKERRQ(ierr);
129*4d9407bcSMatthew G. Knepley   }
130*4d9407bcSMatthew G. Knepley   if (subdm) {
131*4d9407bcSMatthew G. Knepley     PetscSection subsection;
132*4d9407bcSMatthew G. Knepley     PetscBool    haveNull = PETSC_FALSE;
133*4d9407bcSMatthew G. Knepley     PetscInt     f, nf = 0;
134*4d9407bcSMatthew G. Knepley 
135*4d9407bcSMatthew G. Knepley     ierr = PetscSectionCreateSubsection(section, numFields, fields, &subsection);CHKERRQ(ierr);
136*4d9407bcSMatthew G. Knepley     ierr = DMSetDefaultSection(*subdm, subsection);CHKERRQ(ierr);
137*4d9407bcSMatthew G. Knepley     ierr = PetscSectionDestroy(&subsection);CHKERRQ(ierr);
138*4d9407bcSMatthew G. Knepley     for (f = 0; f < numFields; ++f) {
139*4d9407bcSMatthew G. Knepley       (*subdm)->nullspaceConstructors[f] = dm->nullspaceConstructors[fields[f]];
140*4d9407bcSMatthew G. Knepley       if ((*subdm)->nullspaceConstructors[f]) {
141*4d9407bcSMatthew G. Knepley         haveNull = PETSC_TRUE;
142*4d9407bcSMatthew G. Knepley         nf       = f;
143*4d9407bcSMatthew G. Knepley       }
144*4d9407bcSMatthew G. Knepley     }
145*4d9407bcSMatthew G. Knepley     if (haveNull) {
146*4d9407bcSMatthew G. Knepley       MatNullSpace nullSpace;
147*4d9407bcSMatthew G. Knepley 
148*4d9407bcSMatthew G. Knepley       ierr = (*(*subdm)->nullspaceConstructors[nf])(*subdm, nf, &nullSpace);CHKERRQ(ierr);
149*4d9407bcSMatthew G. Knepley       ierr = PetscObjectCompose((PetscObject) *is, "nullspace", (PetscObject) nullSpace);CHKERRQ(ierr);
150*4d9407bcSMatthew G. Knepley       ierr = MatNullSpaceDestroy(&nullSpace);CHKERRQ(ierr);
151*4d9407bcSMatthew G. Knepley     }
152*4d9407bcSMatthew G. Knepley     if (dm->fields) {
153*4d9407bcSMatthew G. Knepley       if (nF != dm->numFields) SETERRQ2(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "The number of DM fields %d does not match the number of Section fields %d", dm->numFields, nF);
154*4d9407bcSMatthew G. Knepley       ierr = DMSetNumFields(*subdm, numFields);CHKERRQ(ierr);
155*4d9407bcSMatthew G. Knepley       for (f = 0; f < numFields; ++f) {
156*4d9407bcSMatthew G. Knepley         ierr = PetscObjectListDuplicate(dm->fields[fields[f]]->olist, &(*subdm)->fields[f]->olist);CHKERRQ(ierr);
157*4d9407bcSMatthew G. Knepley       }
158*4d9407bcSMatthew G. Knepley       if (numFields == 1) {
159*4d9407bcSMatthew G. Knepley         MatNullSpace space;
160*4d9407bcSMatthew G. Knepley         Mat          pmat;
161*4d9407bcSMatthew G. Knepley 
162*4d9407bcSMatthew G. Knepley         ierr = PetscObjectQuery((*subdm)->fields[0], "nullspace", (PetscObject*) &space);CHKERRQ(ierr);
163*4d9407bcSMatthew G. Knepley         if (space) {ierr = PetscObjectCompose((PetscObject) *is, "nullspace", (PetscObject) space);CHKERRQ(ierr);}
164*4d9407bcSMatthew G. Knepley         ierr = PetscObjectQuery((*subdm)->fields[0], "nearnullspace", (PetscObject*) &space);CHKERRQ(ierr);
165*4d9407bcSMatthew G. Knepley         if (space) {ierr = PetscObjectCompose((PetscObject) *is, "nearnullspace", (PetscObject) space);CHKERRQ(ierr);}
166*4d9407bcSMatthew G. Knepley         ierr = PetscObjectQuery((*subdm)->fields[0], "pmat", (PetscObject*) &pmat);CHKERRQ(ierr);
167*4d9407bcSMatthew G. Knepley         if (pmat) {ierr = PetscObjectCompose((PetscObject) *is, "pmat", (PetscObject) pmat);CHKERRQ(ierr);}
168*4d9407bcSMatthew G. Knepley       }
169*4d9407bcSMatthew G. Knepley     }
170*4d9407bcSMatthew G. Knepley   }
171*4d9407bcSMatthew G. Knepley   PetscFunctionReturn(0);
172*4d9407bcSMatthew G. Knepley }
173