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