xref: /petsc/src/dm/interface/dmi.c (revision fad22124558c3e9d4037150e6d1a55fbc78e6dda)
111689aebSJed Brown #include <petsc-private/dmimpl.h>     /*I      "petscdm.h"     I*/
211689aebSJed Brown 
311689aebSJed Brown extern PetscErrorCode VecView_Seq(Vec, PetscViewer);
411689aebSJed Brown extern PetscErrorCode VecView_MPI(Vec, PetscViewer);
511689aebSJed Brown 
611689aebSJed Brown #undef __FUNCT__
711689aebSJed Brown #define __FUNCT__ "DMCreateGlobalVector_Section_Private"
811689aebSJed Brown PetscErrorCode DMCreateGlobalVector_Section_Private(DM dm,Vec *vec)
911689aebSJed Brown {
1011689aebSJed Brown   PetscSection   gSection;
1111689aebSJed Brown   PetscInt       localSize, blockSize = -1, pStart, pEnd, p;
12*fad22124SMatthew G Knepley   PetscErrorCode ierr;
1311689aebSJed Brown 
1411689aebSJed Brown   PetscFunctionBegin;
1511689aebSJed Brown   ierr = DMGetDefaultGlobalSection(dm, &gSection);CHKERRQ(ierr);
16*fad22124SMatthew G Knepley   ierr = PetscSectionGetChart(gSection, &pStart, &pEnd);CHKERRQ(ierr);
1711689aebSJed Brown   for (p = pStart; p < pEnd; ++p) {
1811689aebSJed Brown     PetscInt dof, cdof;
1911689aebSJed Brown 
20*fad22124SMatthew G Knepley     ierr = PetscSectionGetDof(gSection, p, &dof);CHKERRQ(ierr);
21*fad22124SMatthew G Knepley     ierr = PetscSectionGetConstraintDof(gSection, p, &cdof);CHKERRQ(ierr);
2211689aebSJed Brown     if ((blockSize < 0) && (dof > 0)) blockSize = dof-cdof;
2311689aebSJed Brown     if ((dof > 0) && (dof-cdof != blockSize)) {
2411689aebSJed Brown       blockSize = 1;
2511689aebSJed Brown       break;
2611689aebSJed Brown     }
2711689aebSJed Brown   }
28*fad22124SMatthew G Knepley   ierr = PetscSectionGetConstrainedStorageSize(gSection, &localSize);CHKERRQ(ierr);
2911689aebSJed Brown   if (localSize%blockSize) SETERRQ2(((PetscObject) dm)->comm, PETSC_ERR_ARG_WRONG, "Mismatch between blocksize %d and local storage size %d", blockSize, localSize);
3011689aebSJed Brown   ierr = VecCreate(((PetscObject) dm)->comm, vec);CHKERRQ(ierr);
3111689aebSJed Brown   ierr = VecSetSizes(*vec, localSize, PETSC_DETERMINE);CHKERRQ(ierr);
3211689aebSJed Brown   ierr = VecSetBlockSize(*vec, blockSize);CHKERRQ(ierr);
3311689aebSJed Brown   /* ierr = VecSetType(*vec, dm->vectype);CHKERRQ(ierr); */
3411689aebSJed Brown   ierr = VecSetFromOptions(*vec);CHKERRQ(ierr);
3511689aebSJed Brown   ierr = VecSetDM(*vec, dm);CHKERRQ(ierr);
3611689aebSJed Brown   /* ierr = VecSetLocalToGlobalMapping(*vec, dm->ltogmap);CHKERRQ(ierr); */
3711689aebSJed Brown   /* ierr = VecSetLocalToGlobalMappingBlock(*vec, dm->ltogmapb);CHKERRQ(ierr); */
3811689aebSJed Brown   PetscFunctionReturn(0);
3911689aebSJed Brown }
4011689aebSJed Brown 
4111689aebSJed Brown #undef __FUNCT__
4211689aebSJed Brown #define __FUNCT__ "DMCreateLocalVector_Section_Private"
4311689aebSJed Brown PetscErrorCode DMCreateLocalVector_Section_Private(DM dm,Vec *vec)
4411689aebSJed Brown {
45*fad22124SMatthew G Knepley   PetscSection   section;
4611689aebSJed Brown   PetscInt       localSize, blockSize = -1, pStart, pEnd, p;
47*fad22124SMatthew G Knepley   PetscErrorCode ierr;
4811689aebSJed Brown 
4911689aebSJed Brown   PetscFunctionBegin;
50*fad22124SMatthew G Knepley   ierr = DMGetDefaultSection(dm, &section);CHKERRQ(ierr);
51*fad22124SMatthew G Knepley   ierr = PetscSectionGetChart(section, &pStart, &pEnd);CHKERRQ(ierr);
5211689aebSJed Brown   for (p = pStart; p < pEnd; ++p) {
5311689aebSJed Brown     PetscInt dof;
5411689aebSJed Brown 
55*fad22124SMatthew G Knepley     ierr = PetscSectionGetDof(section, p, &dof);CHKERRQ(ierr);
5611689aebSJed Brown     if ((blockSize < 0) && (dof > 0)) blockSize = dof;
5711689aebSJed Brown     if ((dof > 0) && (dof != blockSize)) {
5811689aebSJed Brown       blockSize = 1;
5911689aebSJed Brown       break;
6011689aebSJed Brown     }
6111689aebSJed Brown   }
62*fad22124SMatthew G Knepley   ierr = PetscSectionGetStorageSize(section, &localSize);CHKERRQ(ierr);
6311689aebSJed Brown   ierr = VecCreate(PETSC_COMM_SELF, vec);CHKERRQ(ierr);
6411689aebSJed Brown   ierr = VecSetSizes(*vec, localSize, localSize);CHKERRQ(ierr);
6511689aebSJed Brown   ierr = VecSetBlockSize(*vec, blockSize);CHKERRQ(ierr);
6611689aebSJed Brown   ierr = VecSetFromOptions(*vec);CHKERRQ(ierr);
6711689aebSJed Brown   ierr = VecSetDM(*vec, dm);CHKERRQ(ierr);
6811689aebSJed Brown   PetscFunctionReturn(0);
6911689aebSJed Brown }
70