xref: /petsc/src/dm/tests/ex10.c (revision 8ebe3e4e9e00d86ece2e9fcd0cc84910b0ad437c)
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,ecompkey,vcompkey;
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   /* Register zero size componets to get compkeys to be used by DMNetworkAddComponent() */
30   ierr = DMNetworkRegisterComponent(network,"ecomp",0,&ecompkey);CHKERRQ(ierr);
31   ierr = DMNetworkRegisterComponent(network,"vcomp",0,&vcompkey);CHKERRQ(ierr);
32 
33   Ne = 2;
34   Ni = 1;
35   nodeOffset = (Ne+Ni)*rank;   /* The global node index of the first node defined on this process */
36 
37   /* There are three nodes on each rank and two edges. The edges only connect nodes on the given rank */
38   nedge = k * Ni;
39 
40   ierr = PetscCalloc1(2*nedge,&edgelist);CHKERRQ(ierr);
41   edgelist[0] = nodeOffset + 0;
42   edgelist[1] = nodeOffset + 2;
43   edgelist[2] = nodeOffset + 1;
44   edgelist[3] = nodeOffset + 2;
45 
46   ierr = DMNetworkSetNumSubNetworks(network,PETSC_DECIDE,1);CHKERRQ(ierr);
47   ierr = DMNetworkAddSubnetwork(network,"Subnetwork 1",nedge,edgelist,NULL);CHKERRQ(ierr);
48   ierr = DMNetworkLayoutSetUp(network);CHKERRQ(ierr);
49 
50   /* Add components and variables for the network */
51   ierr = DMNetworkGetSubnetwork(network,0,&nv,&ne,&nodes,&edges);CHKERRQ(ierr);
52   for (e = 0; e < ne; e++) {
53     /* The edges have no degrees of freedom */
54     ierr = DMNetworkAddComponent(network,edges[e],ecompkey,NULL,1);CHKERRQ(ierr);
55   }
56   for (v = 0; v < nv; v++) {
57     ierr = DMNetworkAddComponent(network,nodes[v],vcompkey,NULL,2);CHKERRQ(ierr);
58   }
59 
60   ierr = DMSetUp(network);CHKERRQ(ierr);
61   ierr = PetscPrintf(PETSC_COMM_WORLD,"Network after DMSetUp:\n");CHKERRQ(ierr);
62   ierr = DMView(network,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
63   ierr = DMNetworkGetPlex(network,&plex);CHKERRQ(ierr);
64   /* ierr = DMView(plex,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); */
65   ierr = DMGetLocalSection(plex,&section);CHKERRQ(ierr);
66   ierr = PetscSectionView(section,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
67 
68   ierr = PetscFree(edgelist);CHKERRQ(ierr);
69 
70   ierr = DMNetworkDistribute(&network,0);CHKERRQ(ierr);
71   ierr = PetscPrintf(PETSC_COMM_WORLD,"Network after DMNetworkDistribute:\n");CHKERRQ(ierr);
72   ierr = DMView(network,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
73   ierr = DMNetworkGetPlex(network,&plex);CHKERRQ(ierr);
74   /* ierr = DMView(plex,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); */
75   ierr = DMGetLocalSection(plex,&section);CHKERRQ(ierr);
76   ierr = PetscSectionView(section,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
77 
78   ierr = DMDestroy(&network);CHKERRQ(ierr);
79   ierr = PetscFinalize();
80   return ierr;
81 }
82 
83 /*TEST
84 
85    build:
86       requires: !complex double
87 
88    test:
89       nsize: 2
90 
91 TEST*/
92