[petsc-users] coordinates of vertices of a cell in 3D DMPlex

Mark Adams mfadams at lbl.gov
Sat Jul 25 08:08:33 CDT 2020


Yea, I did not get all the code you need. Here is an example of making
crddm. I'm not sure if this is all best practices (Matt?)

  /* create coordinate DM */
  ierr = DMClone(dm, &crddm);CHKERRV(ierr);
  ierr = PetscFECreateDefault(PETSC_COMM_SELF, dim, dim, PETSC_FALSE, "",
PETSC_DECIDE, &fe);CHKERRV(ierr);
  // ierr = PetscFECreateLagrange(PETSC_COMM_SELF, dim, dim, PETSC_FALSE,
"", PETSC_DECIDE, &fe);CHKERRQ(ierr);
  ierr = PetscFESetFromOptions(fe);CHKERRV(ierr);
  ierr = DMSetField(crddm, field, NULL, (PetscObject)fe);CHKERRV(ierr);
  ierr = DMCreateDS(crddm);CHKERRV(ierr);
  ierr = PetscFEDestroy(&fe);CHKERRV(ierr);

On Sat, Jul 25, 2020 at 7:40 AM Stefano Zampini <stefano.zampini at gmail.com>
wrote:

> Mark
>
> This will only work if you have a vector space for the function
>
>
> On Jul 25, 2020, at 1:13 PM, Mark Adams <mfadams at lbl.gov> wrote:
>
> I get all the coordinates with this method:
>
> static PetscErrorCode crd_func(PetscInt dim, PetscReal time, const
> PetscReal x[], PetscInt Nf_dummy, PetscScalar *u, void *actx)
> {
>   int i;
>   PetscFunctionBeginUser;
>   for (i = 0; i < dim; ++i) u[i] = x[i];
>   PetscFunctionReturn(0);
> }
>
>   PetscErrorCode (*initu[1])(PetscInt, PetscReal, const PetscReal [],
> PetscInt, PetscScalar [], void *);
>   /* project coordinates to vertices */
>   ierr = DMCreateGlobalVector(crddm, &crd_vec);CHKERRV(ierr);
>   initu[0] = crd_func;
>   ierr = DMProjectFunction(crddm, 0.0, initu, NULL, INSERT_ALL_VALUES,
> crd_vec);CHKERRV(ierr);
>   ierr = VecViewFromOptions(crd_vec, NULL, "-coord_view");CHKERRV(ierr);
>   /* iterate over mesh data and get indices */
>   ierr = VecGetArrayRead(crd_vec,&xx);CHKERRV(ierr);
>   ierr = VecGetLocalSize(rho,&N);CHKERRV(ierr);
>   /* access grid data here */
>   for (p=0;p<N;p++) {
>     for (i=0;i<dim;i++) ..... = xx[p*dim+i];
>     PetscPrintf(PETSC_COMM_SELF,"xx = (%g, %g)\n", xx[p*dim+0],
> xx[p*dim+1]);
>   }
>   ierr = VecRestoreArrayRead(crd_vec,&xx);CHKERRV(ierr);
>   ierr = VecDestroy(&crd_vec);CHKERRV(ierr);
>
> On Sat, Jul 25, 2020 at 4:10 AM Swarnava Ghosh <swarnava89 at gmail.com>
> wrote:
>
>> Dear Petsc users,
>>
>> I had a trivial question about DMPlex. Suppose I have a 3D mesh of
>> tetrahedrons. I want to find out the 3D coordinates of the vertices of a
>> particular cell. What would be the function to do this?
>>
>> Thank you,
>> SG
>>
>
>
-------------- next part --------------
An HTML attachment was scrubbed...
URL: <http://lists.mcs.anl.gov/pipermail/petsc-users/attachments/20200725/076cf90d/attachment-0001.html>


More information about the petsc-users mailing list