[petsc-users] Gauss-Lobatto-Legendre Element Gradient -- Caught signal number 11 SEGV

Matthew Knepley knepley at gmail.com
Thu Jun 1 06:24:02 CDT 2023


On Thu, Jun 1, 2023 at 1:46 AM Duan Junming via petsc-users <
petsc-users at mcs.anl.gov> wrote:

> Dear all,
>
>
> I have a simple demo code attached below, which gives a segmentation
> violation error.
>
> Can you help me with this problem? I think the problem is due to the
> destroy function.
>
> I am using version 3.19.2 with debugging.
>

Yes, the check is wrong. I will fix it. For now, you can just pass in the
transpose argument as well.

  Thanks,

    Matt


> #include <petsc.h>
>
> static char help[] = "test.\n";
> int main(int argc, char *argv[]) {
>   PetscCall(PetscInitialize(&argc, &argv, 0, help));
>   PetscScalar *nodes;
>   PetscScalar *weights;
>   PetscScalar **diff;
>   PetscInt n = 3;
>   PetscCall(PetscMalloc2(n, &nodes, n, &weights));
>   PetscCall(PetscDTGaussLobattoLegendreQuadrature(n,
> PETSCGAUSSLOBATTOLEGENDRE_VIA_LINEAR_ALGEBRA, nodes, weights));
>   PetscCall(PetscGaussLobattoLegendreElementGradientCreate(n, nodes,
> weights, &diff, NULL));
>   PetscCall(PetscGaussLobattoLegendreElementGradientDestroy(n, nodes,
> weights, &diff, NULL));
>   PetscCall(PetscFree2(nodes, weights));
>   PetscCall(PetscFinalize());
>   return 0;
> }
>
>
> Junming
>
>

-- 
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://www.cse.buffalo.edu/~knepley/ <http://www.cse.buffalo.edu/~knepley/>
-------------- next part --------------
An HTML attachment was scrubbed...
URL: <http://lists.mcs.anl.gov/pipermail/petsc-users/attachments/20230601/24bc801d/attachment.html>


More information about the petsc-users mailing list