On Fri, Feb 23, 2024 at 9:54 AM 袁煕 <[email protected]> wrote:
> Hello, I cannot fetch right values of a local vector after a second time
> calling DMSetLocalSection. My program runs like ==================== call
> PetscSectionSetup(section1,ierr) call
> DMSetLocalSection(dm_mesh,section1,ierr)call
> PetscSectionDestroy(section1,ierr)call
>
> ZjQcmQRYFpfptBannerStart
> This Message Is From an External Sender
> This message came from outside your organization.
>
> ZjQcmQRYFpfptBannerEnd
> Hello,
>
> I cannot fetch right values of a local vector after a second time calling
> DMSetLocalSection. My program runs like
>
> ====================
> call PetscSectionSetup(section1,ierr)
> call DMSetLocalSection(dm_mesh,section1,ierr)
> call PetscSectionDestroy(section1,ierr)
>
> call DMCreateGlobalVector(dm_mesh,vec_a,ierr)
>
> call VecDestroy(vec_a,ierr)
> call DMClearLocalVectors(dm_mesh, ierr);
>
You also need to wipe out the SF,
call DMSetSectionSF(dm, PETSC_NULL_PETSCSF, ierr)
> call PetscSectionSetup(section2,ierr)
> call DMSetLocalSection(dm_mesh,section2,ierr)
>
> call DMCreateGlobalVector(dm_mesh,vec_a,ierr)
> call DMGetLocalVector(dm_mesh,lvec,ierr);
> call DMGlobalToLocal(dm_mesh,vec_a,INSERT_VALUES,lvec,ierr)
> call VecGetArrayReadF90(lvec,tempv,ierr)
> call VecGetArrayReadF90(vec_a,tempv1,ierr)
> ========================
>
> In this program, I got two Fortran style arrays, *tempv* and *tempv1*,
> from a Vec *vec_a*. When only 1 cpu is used, those two arrays should be
> completely the same. Am I right? But I found there are differences as
> follows:
>
> * In array *tempv*, dofs newly introduced by section2 are zero valued.
> Looks like those values are ignored by DMGlobalToLocal. But in array
> *tempv1*, all dofs values are set correctly.
>
> Did I do something wrong in the above code?
>
We do not usually reset the section on a DM. Instead, when we want to use a
different function space,
we clone the DM, using DMClone(), and make a new function space on the new
DM.
Thanks,
MAtt
Much thanks for your help.
>
> X. Yuan, Ph.D. in Solid Mechanics
>
--
What most experimenters take for granted before they begin their
experiments is infinitely more interesting than any results to which their
experiments lead.
-- Norbert Wiener
https://urldefense.us/v3/__https://www.cse.buffalo.edu/*knepley/__;fg!!G_uCfscf7eWS!bAsyX_068Z0c_UIQmCwlZXop22TY8QcA8rRgki-yA9aSOze_tjnsWkrWUfzE8KpRNzzHsqLXQwFSAMO0WIRk$
<https://urldefense.us/v3/__http://www.cse.buffalo.edu/*knepley/__;fg!!G_uCfscf7eWS!bAsyX_068Z0c_UIQmCwlZXop22TY8QcA8rRgki-yA9aSOze_tjnsWkrWUfzE8KpRNzzHsqLXQwFSADY8jt9N$
>