[petsc-users] distribute and cells mapping.

Olivier Bonnefon olivier.bonnefon at avignon.inra.fr
Wed Aug 28 10:01:53 CDT 2013


On 08/28/2013 10:28 AM, Olivier Bonnefon wrote:
> On 08/23/2013 04:42 PM, Matthew Knepley wrote:
>> On Fri, Aug 23, 2013 at 9:35 AM, Olivier Bonnefon 
>> <olivier.bonnefon at avignon.inra.fr 
>> <mailto:olivier.bonnefon at avignon.inra.fr>> wrote:
>>
>>     Hello,
>>
>>     Thanks for your answers, I'm now able to import and distribute a
>>     mesh:
>>
>>
>> You might simplify this to
>>
>> if (rank) {obNbCells = 0; obNbVertex = 0;}
>> ierr = 
>> DMPlexCreateFromCellList(comm,dim,obNbCells,obNbVertex,3,0,obCells,2,obVertex,dm);CHKERRQ(ierr);
>>
>>     if (!rank){
>>             ierr =
>>     DMPlexCreateFromCellList(comm,dim,obNbCells,obNbVertex,3,0,obCells,2,obVertex,dm);CHKERRQ(ierr);
>>              for (i=0;i<obNbBound;i++){
>>                ierr =DMPlexSetLabelValue(*dm, "marker",
>>     obBoundary[i]+obNbCells, 1);CHKERRQ(ierr);
>>              }
>>     }else {
>>               ierr =
>>     DMPlexCreateFromCellList(comm,dim,0,0,3,0,obCells,2,obVertex,dm);CHKERRQ(ierr);
>>     }
>>
>>     ierr = DMPlexDistribute(*dm, partitioner, 0,
>>     &distributedMesh);CHKERRQ(ierr);
>>     if (distributedMesh) {
>>           ierr = DMDestroy(dm);CHKERRQ(ierr);
>>           *dm  = distributedMesh;
>>     }
>>
>>     Is it possible to known the resulting partition ? ie, What is the
>>     mapping between the initial cell number and the local cell (used
>>     in DMPlexComputeResidualFEM)?
>>     I need this to write an efficient implementation of the FEM
>>     struct functions f0 and g0, space depending.
>>
>>
>> Yes, but I really do not think you want to do things that way. I am 
>> assuming you want different material models or something
>> in different places. The way I envision that is using a DMLabel to 
>> mark up parts of the domain. All labels are automatically
>> distributed with the mesh. Is that what you want?
> Hello,
>
> It is exactly what I need: I'm mobilized a landscape, and the 
> parameters of the model depend of the type of crop. Therefore, I have 
> created a label for each type of crop and I have labeled each triangle 
> with the corresponding label:
>
>  for (i=0;i<obNbCells;i++){
>      if (labelCells[i]==1){
>         ierr =DMPlexSetLabelValue(*dm, "marker1", i, 1);CHKERRQ(ierr);
>      }else{
>         ierr =DMPlexSetLabelValue(*dm, "marker2", i, 1);CHKERRQ(ierr);
>      }
>  }
>
> So, I'm able to mark the triangles, but I'm not able to get this label 
> in the plugin "fem.f0Funcs" and "fem.g0Funcs": These plugins are 
> called by looping on the triangles in the function 
> "FEMIntegrateResidualBatch", but the dm is not available, so I can't 
> use the functions DMPlexGetLabel, DMLabelGetStratumSize and 
> DMLabelGetStratumIS. What is the good way to get the labels in the 
> user plugins of the fem struct ?
>
>
> Thanks a lot for your help.
>
> Olivier B
Hello,

This is the solution I implemented to get the label level in the plugins 
"fem.f0Funcs" and "fem.g0Funcs":

  I need the DM and the index element, so i do:
1) I  add some static variables:
static DM * spDM[128];
static int scurElem[128];
2) I overload the DMPlexComputeJacobianFEM with :
    PetscErrorCode MY_DMPlexComputeJacobianFEM(DM dm, Vec X, Mat Jac, 
Mat JacP, MatStructure *str,void *user)
{

  PetscMPIInt rank;
  PetscErrorCode ierr;

  PetscFunctionBeginUser;
  ierr = MPI_Comm_rank(PETSC_COMM_WORLD, &rank);CHKERRQ(ierr);
  spDM[rank]=&dm;
  PetscFunctionReturn(DMPlexComputeJacobianFEM(dm, X, Jac,JacP, str,user));

}
3) overload FEMIntegrateResidualBatch adding code:
.
.
   for (e = 0; e < Ne; ++e) {
     scurElem[rank]=e;//added ligne
.
.

So that, I can get the label level using  DMPlexHasLabel and DMLabelGetValue

I'm sure this solution is awful, and works only in this version, but i 
didn't find a better way to get the label in the plugins fem struc. Do 
you know the correct way to do that ??

Thanks,
Olivier B
>>
>>   Thanks,
>>
>>      Matt
>>
>>     Regards,
>>
>>     Olivier B
>>
>>     -- 
>>     Olivier Bonnefon
>>     INRA PACA-Avignon, Unité BioSP
>>     Tel: +33 (0)4 32 72 21 58 <tel:%2B33%20%280%294%2032%2072%2021%2058>
>>
>>
>>
>>
>> -- 
>> 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
>
>
> -- 
> Olivier Bonnefon
> INRA PACA-Avignon, Unité BioSP
> Tel: +33 (0)4 32 72 21 58


-- 
Olivier Bonnefon
INRA PACA-Avignon, Unité BioSP
Tel: +33 (0)4 32 72 21 58

-------------- next part --------------
An HTML attachment was scrubbed...
URL: <http://lists.mcs.anl.gov/pipermail/petsc-users/attachments/20130828/bc0d56ee/attachment.html>


More information about the petsc-users mailing list