[petsc-users] dmda and structures

Gideon Simpson gideon.simpson at icloud.com
Tue Dec 31 20:06:04 CST 2013


I believe I've gotten things to work in both cases, but I wanted to clarify what the substantial difference between DMDAVecGetArray and DMDAVecGetArrayDOF.  Is it just that it allows me to do indexing of the form realization_array[i][j] instead of  realization_array[i * num_dof +j], or is there something more here, in terms of how things are distributed?

-gideon

On Dec 29, 2013, at 2:53 PM, Barry Smith <bsmith at mcs.anl.gov> wrote:

> 
> On Dec 29, 2013, at 1:38 PM, Gideon Simpson <gideon.simpson at icloud.com> wrote:
> 
>> How would that work?
> 
>   I’ve modified your code fragment below with how I think you could do it.
> 
>> I took Jed's advice, and set it up so that I have access via:
>> 
>> realization_array[ realization_dof * i + 1 +j] = jth point of ith realization
>> 
>> and then I further used a pointer:
>> 
>> PetscScalar *x;
>> 
>> x = & realization_array[ realization_dof * i + 1];
>> 
>> so that I have "clean" access to the elements of each realization through x[j].
>> 
>> -gideon
>> 
>> On Dec 29, 2013, at 2:33 PM, Barry Smith <bsmith at mcs.anl.gov> wrote:
>> 
>>> 
>>> DMDAVecGetArrayDOF() might serve your needs. You would access with realization_array[i][j]  for the ith grid point and j th realization
>>> 
>>> Barry
>>> 
>>> 
>>> 
>>> On Dec 28, 2013, at 7:44 PM, Gideon Simpson <gideon.simpson at icloud.com> wrote:
>>> 
>>>> I'm trying to use 1D DMDA's to manage distributed monte carlo computations, and I've got the following problem.  If I know, in advance, how many floating points are needed for each realization, everything is fine.  But, if I want to be able to set this at run time, I can't get it to work.  What I have is:
>>>> 
>>>> 
>>>> then, in the main part of the code:
>>>> 
>>>> Realization *realization_array;
>>>> DM sim_data_da;
>>>> Vec sim_data_vec;
>>>> 
>>>> PetscInt xsize=10, realization_dof, i;
>>>> PetscOptionsGetInt(NULL, "-xsize", &xsize, NULL);
>>>> 
>>>> realization_dof	  = xsize + 1;
>>>> 
>>>> 
>>>> DMDACreate1d(PETSC_COMM_WORLD,DMDA_BOUNDARY_NONE, batch_size, realization_dof, 0, NULL, &sim_data_da);
>>>> DMCreateGlobalVector(sim_data_da,&sim_data_vec);
>>>> 
> PetscScalar **realization_array;
> 
>>>> DMDAVecGetArrayDOF(sim_data_da, sim_data_vec, & realization_array);
> 
>    something = realization_array[i][j];     here i varies from xs to xs+xm-1  and j varies from 0 to realization_dof-1
> 
>>>> 
>>>> 
>>>> -gideon
>>>> 
>>> 
>> 
> 



More information about the petsc-users mailing list