xref: /petsc/src/dm/tests/ex10.c (revision 8dbb0df66754dee4fb72dff2ad56e76db1e6b7c7)
1 /*
2     Simple example demonstrating that creating a one sub-network DMNetwork in parallel.
3 */
4 
5 #include <petscdmnetwork.h>
6 
7 int main(int argc,char ** argv)
8 {
9   PetscErrorCode    ierr;
10   DM                network;
11   PetscMPIInt       size,rank;
12   MPI_Comm          comm;
13   PetscInt          e,ne,nv,v;
14   PetscInt          *edgelist = NULL;
15   const PetscInt    *nodes,*edges;
16   DM                plex;
17   PetscSection      section;
18   PetscInt          Ne,Ni;
19   PetscInt          nodeOffset,k = 2,nedge;
20 
21   ierr = PetscInitialize(&argc,&argv,NULL,NULL);if (ierr) return ierr;
22   ierr = PetscOptionsSetValue(NULL,"-petscpartitioner_use_vertex_weights","No");CHKERRQ(ierr);
23   comm = PETSC_COMM_WORLD;
24   ierr = MPI_Comm_rank(comm,&rank);CHKERRMPI(ierr);
25   ierr = MPI_Comm_size(comm,&size);CHKERRMPI(ierr);
26 
27   ierr = DMNetworkCreate(PETSC_COMM_WORLD,&network);CHKERRQ(ierr);
28 
29   Ne = 2;
30   Ni = 1;
31   nodeOffset = (Ne+Ni)*rank;   /* The global node index of the first node defined on this process */
32 
33   /* There are three nodes on each rank and two edges. The edges only connect nodes on the given rank */
34   nedge = k * Ni;
35 
36   ierr = PetscCalloc1(2*nedge,&edgelist);CHKERRQ(ierr);
37   edgelist[0] = nodeOffset + 0;
38   edgelist[1] = nodeOffset + 2;
39   edgelist[2] = nodeOffset + 1;
40   edgelist[3] = nodeOffset + 2;
41 
42   ierr = DMNetworkSetNumSubNetworks(network,PETSC_DECIDE,1);CHKERRQ(ierr);
43   ierr = DMNetworkAddSubnetwork(network,"Subnetwork 1",nedge,edgelist,NULL);CHKERRQ(ierr);
44   ierr = DMNetworkLayoutSetUp(network);CHKERRQ(ierr);
45 
46   /* Add components and variables for the network */
47   ierr = DMNetworkGetSubnetwork(network,0,&nv,&ne,&nodes,&edges);CHKERRQ(ierr);
48   for (e = 0; e < ne; e++) {
49     /* The edges have no degrees of freedom */
50     ierr = DMNetworkAddComponent(network,edges[e],-1,NULL,1);CHKERRQ(ierr);
51   }
52   for (v = 0; v < nv; v++) {
53     ierr = DMNetworkAddComponent(network,nodes[v],-1,NULL,2);CHKERRQ(ierr);
54   }
55 
56   ierr = DMSetUp(network);CHKERRQ(ierr);
57   ierr = PetscPrintf(PETSC_COMM_WORLD,"Network after DMSetUp:\n");CHKERRQ(ierr);
58   ierr = DMView(network,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
59   ierr = DMNetworkGetPlex(network,&plex);CHKERRQ(ierr);
60   /* ierr = DMView(plex,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); */
61   ierr = DMGetLocalSection(plex,&section);CHKERRQ(ierr);
62   ierr = PetscSectionView(section,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
63 
64   ierr = PetscFree(edgelist);CHKERRQ(ierr);
65 
66   ierr = DMNetworkDistribute(&network,0);CHKERRQ(ierr);
67   ierr = PetscPrintf(PETSC_COMM_WORLD,"Network after DMNetworkDistribute:\n");CHKERRQ(ierr);
68   ierr = DMView(network,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
69   ierr = DMNetworkGetPlex(network,&plex);CHKERRQ(ierr);
70   /* ierr = DMView(plex,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); */
71   ierr = DMGetLocalSection(plex,&section);CHKERRQ(ierr);
72   ierr = PetscSectionView(section,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
73 
74   ierr = DMDestroy(&network);CHKERRQ(ierr);
75   ierr = PetscFinalize();
76   return ierr;
77 }
78 
79 /*TEST
80 
81    build:
82       requires: !complex double
83 
84    test:
85       nsize: 2
86 
87 TEST*/
88