[petsc-users] A question of DMPlexGetAdjacency
Matthew Knepley
knepley at gmail.com
Tue Nov 15 06:09:45 CST 2016
On Tue, Nov 15, 2016 at 6:07 AM, leejearl <leejearl at 126.com> wrote:
> Hi, Matt:
> Thank you for your kind reply. Since it is a low level routine, is
> there any method which I can used for creating a
> dynamic array?
>
I am not sure what you are asking for?
Thanks,
Matt
> Thanks
> leejearl
>
> On 2016年11月15日 16:20, Matthew Knepley wrote:
>
> On Tue, Nov 15, 2016 at 1:00 AM, leejearl <leejearl at 126.com> wrote:
>
>> Hi Lawrence:
>>
>> Thank you for your prompt reply. It works well, but I encounter
>> another problem.
>>
>> The code is as follows:
>>
>> PetscScalar *arrayConVar, *arrayDeltaX, *arrayDeltaY;
>>
>> ierr = DMPlexGetAdjacency(dm, c, &adjSize, &adj); CHKERRQ(ierr);
>>
>> /*PetscSynchronizedPrintf(PETSC_COMM_WORLD, "adjSize = %d\n",
>> adjSize);*/
>> adjSizeInterior = -1;
>> for(i = 0; i < adjSize; i++)
>> {
>> if(adj[i] < cEndInterior)
>> {
>> adjSizeInterior++;
>> /*PetscSynchronizedPrintf(PETSC_COMM_WORLD, "adj[%d]=
>> %d\n", i, adj[i]);*/
>> }
>> }
>> /*
>> PetscSynchronizedPrintf(PETSC_COMM_WORLD, "adjSizeInterior =
>> %d\n", adjSizeInterior);
>> PetscSynchronizedFlush(PETSC_COMM_WORLD,PETSC_STDOUT);
>> */
>>
>> ierr = PetscMalloc3(adjSizeInterior*4, &arrayConVar,
>> adjSizeInterior, &arrayDeltaX, adjSizeInterior, &arrayDeltaY);
>> CHKERRQ(ierr);
>>
>> PetscFree3(arrayConVar, arrayDeltaX, arrayDeltaY); CHKERRQ(ierr);
>> ierr = PetscFree(adj); CHKERRQ(ierr);
>>
>> The error messages are like as :
>>
>> [0]PETSC ERROR: --------------------- Error Message
>> --------------------------------------------------------------
>> [0]PETSC ERROR: Petsc has generated inconsistent data
>> [0]PETSC ERROR: Invalid mesh exceeded adjacency allocation (0)
>>
>
> You are resetting the size to 0 in your loop. It needs to be set to
> PETSC_DETERMINE.
>
> This is a very low level routine with complicated memory management for
> efficiency. It is not intended for
> most user code. However, here is an example:
>
> https://bitbucket.org/petsc/petsc/src/f98bafe09724ec53686c21e51e969e
> 41f10d977b/src/dm/impls/plex/plexdistribute.c?at=master&
> fileviewer=file-view-default#plexdistribute.c-570
>
> Thanks,
>
> Matt
>
>
>> [0]PETSC ERROR: See http://www.mcs.anl.gov/petsc/documentation/faq.html
>> for trouble shooting.
>> [0]PETSC ERROR: Petsc Development GIT revision: v3.7.4-1864-ga10654c GIT
>> Date: 2016-11-01 22:25:06 -0500
>> [0]PETSC ERROR: ./cavity on a arch-linux2-c-debug named leejearl by
>> leejearl Tue Nov 15 14:48:57 2016
>> [0]PETSC ERROR: Configure options --prefix=/home/leejearl/Install/Petsc
>> --with-mpi-dir=/home/leejearl/Install/openmpi/1.8.5_gnu
>> --download-exodusii=../externalpackages/exodus-5.24.tar.bz2
>> --download-netcdf=../externalpackages/netcdf-4.3.2.tar.gz
>> --download-hdf5=../externalpackages/hdf5-1.8.12.tar.gz
>> --download-triangle=../externalpackages/Triangle.tar.gz
>> --download-sowing=../externalpackages/git.sowing.tar.gz
>> [0]PETSC ERROR: #1 DMPlexGetAdjacency_Cone_Internal() line 193 in
>> /home/leejearl/Software/petsc/petsc/src/dm/impls/plex/plexdistribute.c
>> [0]PETSC ERROR: #2 DMPlexGetAdjacency_Internal() line 298 in
>> /home/leejearl/Software/petsc/petsc/src/dm/impls/plex/plexdistribute.c
>> [0]PETSC ERROR: #3 DMPlexGetAdjacency() line 372 in
>> /home/leejearl/Software/petsc/petsc/src/dm/impls/plex/plexdistribute.c
>> [0]PETSC ERROR: #4 CalcConVarGradient() line 126 in
>> /home/leejearl/Desktop/PETSc/GKSCavity_Petsc/cavity.c
>>
>> It is very strange that PetscMalloc3 can arouse such an error. If I
>> comment the statement relative with PetscMalloc3, it is ok! I can not
>> understand what happened.
>>
>> Can anyone figure me out the problems of my codes?
>>
>> I attach the code and the grid, so you can compile in your workstation.
>>
>> Thanks
>>
>> leejearl
>>
>>
>> On 2016年11月14日 23:09, Lawrence Mitchell wrote:
>>
>>> On 14 Nov 2016, at 14:59, leejearl <leejearl at 126.com> wrote:
>>>>
>>>> Hi all:
>>>>
>>>> I am not sure the using of the function "DMPlexGetAdjacency".
>>>>
>>>> My codes are as follows:
>>>>
>>>> PetscInt adj, *adjSize=NULL;
>>>>
>>>> ierr = DMPlexGetAdjacency(dm, p, &adj, &adj);
>>>>
>>>
>>> Your calling sequence is wrong:
>>>
>>> PetscInt adjSize;
>>> PetscInt *adj = NULL;
>>>
>>> ierr = DMPlexGetAdjacency(dm, p, &adjSize, &adj); CHKERRQ(ierr);
>>>
>>> adjSize is now the size of the adj array. You should remember to free
>>> it afterwards:
>>>
>>> /* use adj here *
>>> ...
>>>
>>> ierr = PetscFree(adj); CHKERRQ(ierr);
>>>
>>> Lawrence
>>>
>>
>>
>
>
> --
> 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
>
>
>
--
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
-------------- next part --------------
An HTML attachment was scrubbed...
URL: <http://lists.mcs.anl.gov/pipermail/petsc-users/attachments/20161115/20397288/attachment-0001.html>
More information about the petsc-users
mailing list