:orphan: # DMPlexLocalVectorLoad Loads on-disk vector data into a local vector ## Synopsis ``` #include "petscdmplex.h" PetscErrorCode DMPlexLocalVectorLoad(DM dm, PetscViewer viewer, DM sectiondm, PetscSF sf, Vec vec) ``` Collective ## Input Parameters - ***dm -*** The `DM` that represents the topology - ***viewer -*** The `PetscViewer` that represents the on-disk vector data - ***sectiondm -*** The `DM` that contains the local section on which vec is defined - ***sf -*** The `PetscSF` that migrates the on-disk vector data into vec - ***vec -*** The local vector to set values of ## Notes In general `dm` and `sectiondm` are two different objects, the former carrying the topology and the latter carrying the section, and have been given a topology name and a section name, respectively, with `PetscObjectSetName()`. In practice, however, they can be the same object if it carries both topology and section; in that case the name of the object is used as both the topology name and the section name. ## Typical calling sequence ```none DMCreate(PETSC_COMM_WORLD, &dm); DMSetType(dm, DMPLEX); PetscObjectSetName((PetscObject)dm, "topologydm_name"); DMPlexTopologyLoad(dm, viewer, &sfX); DMClone(dm, §iondm); PetscObjectSetName((PetscObject)sectiondm, "sectiondm_name"); DMPlexSectionLoad(dm, viewer, sectiondm, sfX, NULL, &lsf); DMGetLocalVector(sectiondm, &vec); PetscObjectSetName((PetscObject)vec, "vec_name"); DMPlexLocalVectorLoad(dm, viewer, sectiondm, lsf, vec); DMRestoreLocalVector(sectiondm, &vec); PetscSFDestroy(&lsf); PetscSFDestroy(&sfX); DMDestroy(§iondm); DMDestroy(&dm); ``` ## See Also [](ch_unstructured), `DM`, `DMPLEX`, `DMPlexTopologyLoad()`, `DMPlexSectionLoad()`, `DMPlexGlobalVectorLoad()`, `DMPlexGlobalVectorView()`, `DMPlexLocalVectorView()`, `PetscSF`, `PetscViewer` ## Level advanced ## Location src/dm/impls/plex/plex.c --- [Edit on GitLab](https://gitlab.com/petsc/petsc/-/edit/release/src/dm/impls/plex/plex.c) [Index of all DMPlex routines](index.md) [Table of Contents for all manual pages](/manualpages/index.md) [Index of all manual pages](/manualpages/singleindex.md)