[petsc-users] difference between DMDAVecGetArrayDOF and DMDAVecGetArray?

Matthew Knepley knepley at gmail.com
Sat Feb 28 13:59:22 CST 2015


On Sat, Feb 28, 2015 at 1:56 PM, Gideon Simpson <gideon.simpson at gmail.com>
wrote:

> Wait, then why is my code working when I use the DOF routines?  The number
> of degrees of freedom per DA vertex is constant across the DA, but that
> value is not set until run time.  In other words, what’s wrong with the
> code:
>
> PetscOptionsGetInt(NULL,"-n_samples",&n_samples,NULL);
> PetscOptionsGetInt(NULL,"-n_points",&n_points,NULL);
>
> DMDACreate1d(PETSC_COMM_WORLD,DM_BOUNDARY_NONE, n_samples, n_points, 0,
> NULL, &da);
>

Nothing is wrong with that. I meant that you must have the smae  number on
every vertex.

  Matt


> -gideon
>
> On Feb 28, 2015, at 2:54 PM, Matthew Knepley <knepley at gmail.com> wrote:
>
> On Sat, Feb 28, 2015 at 1:50 PM, Gideon Simpson <gideon.simpson at gmail.com>
> wrote:
>
>> Supposing that I do not have a priori information as to the number of
>> degrees of freedom per DA vertex; I want the number of mesh points along
>> each sample path to be variable.  Hence, I can’t really use a statically
>> defined structure as you suggest.  In that case, are the DOF routines the
>> only option?
>>
>
> DMDA just plain does not support that. It is close, but everything is not
> hooked up right. If you have
> variable numbers of unknowns per site, you can
>
>   1) Fix DMDA with our help (demands programming)
>
>   2) Use DMPlex (demands reading and maybe looking at code)
>
> Thanks,
>
>     Matt
>
>
>> -gideon
>>
>> On Feb 28, 2015, at 2:42 PM, Matthew Knepley <knepley at gmail.com> wrote:
>>
>> On Sat, Feb 28, 2015 at 1:35 PM, Gideon Simpson <gideon.simpson at gmail.com
>> > wrote:
>>
>>> I’m having some trouble understanding what the difference between these
>>> two routines are, though I am finding that there certainly is a
>>> difference.  I have the following monte carlo problem.  I am generating
>>> n_sample paths each of length n_points, and storing them in a 1D DA:
>>>
>>> DMDACreate1d(PETSC_COMM_WORLD,DM_BOUNDARY_NONE, n_samples, n_points, 0,
>>> NULL, &da);
>>> DMCreateGlobalVector(da,&paths_vec);
>>>
>>> When I then go to access them,
>>>
>>> PetscScalar **u_array;
>>>
>>> I find that:
>>>
>>> DMDAVecGetArrayDOF(da, paths_vec, &u_array);
>>>
>>> works as exepected, in that u_array[i] is a pointer to the first index
>>> of the i-th sample path, but if I call:
>>>
>>> DMDAVecGetArray(da, paths_vec, &u_array);
>>>
>>> u_array[i] is something else, and my attempts to manipulate it result in
>>> segmentation faults, even though the code compiles and builds.
>>
>>
>> Suppose that you have 4 PetscScalar values at each vertex of the 1D DMDA.
>> If you use
>>
>>   PetscScalar **u;
>>
>>   DMDAVecGetArrayDOF(da, uVec, &u);
>>
>>   u[i][2] /* refers to 3rd scalar on vertex i */
>>
>> On the other hand you could use
>>
>> typedef struct {
>>   PetscScalar a, b, c, d;
>> } Vals;
>>
>>  Vals *u;
>>
>>  DMDAVecGetArray(da, uVec, &u);
>>
>>  u[i].c /* refers to the same value as above */
>>
>> Basically the DOF version gives you an extra level of indirection for the
>> components.
>>
>>   Thanks,
>>
>>      Matt
>>
>>
>>>
>>> -gideon
>>>
>>>
>>
>>
>> --
>> 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
>
>
>


-- 
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/20150228/8a291bce/attachment.html>


More information about the petsc-users mailing list