1af0996ceSBarry Smith #include <petsc/private/dmnetworkimpl.h> /*I "petscdmnetwork.h" I*/ 25f2c45f1SShri Abhyankar #include <petscdmplex.h> 35f2c45f1SShri Abhyankar #include <petscsf.h> 45f2c45f1SShri Abhyankar 55f2c45f1SShri Abhyankar /*@ 6556ed216SShri Abhyankar DMNetworkGetPlex - Gets the Plex DM associated with this network DM 7556ed216SShri Abhyankar 8556ed216SShri Abhyankar Not collective 9556ed216SShri Abhyankar 10556ed216SShri Abhyankar Input Parameters: 11556ed216SShri Abhyankar + netdm - the dm object 12556ed216SShri Abhyankar - plexmdm - the plex dm object 13556ed216SShri Abhyankar 14556ed216SShri Abhyankar Level: Advanced 15556ed216SShri Abhyankar 16556ed216SShri Abhyankar .seealso: DMNetworkCreate() 17556ed216SShri Abhyankar @*/ 18556ed216SShri Abhyankar PetscErrorCode DMNetworkGetPlex(DM netdm, DM *plexdm) 19556ed216SShri Abhyankar { 20556ed216SShri Abhyankar DM_Network *network = (DM_Network*) netdm->data; 21556ed216SShri Abhyankar 22556ed216SShri Abhyankar PetscFunctionBegin; 23556ed216SShri Abhyankar *plexdm = network->plex; 24556ed216SShri Abhyankar PetscFunctionReturn(0); 25556ed216SShri Abhyankar } 26556ed216SShri Abhyankar 27556ed216SShri Abhyankar /*@ 28e2aaf10cSShri Abhyankar DMNetworkSetSizes - Sets the number of subnetworks,local and global vertices and edges for each subnetwork. 295f2c45f1SShri Abhyankar 305f2c45f1SShri Abhyankar Collective on DM 315f2c45f1SShri Abhyankar 325f2c45f1SShri Abhyankar Input Parameters: 335f2c45f1SShri Abhyankar + dm - the dm object 34e2aaf10cSShri Abhyankar . Nsubnet - number of subnetworks 35e2aaf10cSShri Abhyankar . nV - number of local vertices for each subnetwork 36e2aaf10cSShri Abhyankar . nE - number of local edges for each subnetwork 37e2aaf10cSShri Abhyankar . NV - number of global vertices (or PETSC_DETERMINE) for each subnetwork 38e2aaf10cSShri Abhyankar - NE - number of global edges (or PETSC_DETERMINE) for each subnetwork 395f2c45f1SShri Abhyankar 405f2c45f1SShri Abhyankar Notes 415f2c45f1SShri Abhyankar If one processor calls this with NV (NE) of PETSC_DECIDE then all processors must, otherwise the prgram will hang. 425f2c45f1SShri Abhyankar 435f2c45f1SShri Abhyankar You cannot change the sizes once they have been set 445f2c45f1SShri Abhyankar 451b266c99SBarry Smith Level: intermediate 461b266c99SBarry Smith 471b266c99SBarry Smith .seealso: DMNetworkCreate() 485f2c45f1SShri Abhyankar @*/ 49e2aaf10cSShri Abhyankar PetscErrorCode DMNetworkSetSizes(DM dm, PetscInt Nsubnet, PetscInt nV[], PetscInt nE[], PetscInt NV[], PetscInt NE[]) 505f2c45f1SShri Abhyankar { 515f2c45f1SShri Abhyankar PetscErrorCode ierr; 525f2c45f1SShri Abhyankar DM_Network *network = (DM_Network*) dm->data; 53e2aaf10cSShri Abhyankar PetscInt a[2],b[2],i; 545f2c45f1SShri Abhyankar 555f2c45f1SShri Abhyankar PetscFunctionBegin; 565f2c45f1SShri Abhyankar PetscValidHeaderSpecific(dm,DM_CLASSID,1); 57e2aaf10cSShri Abhyankar if (Nsubnet <= 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Number of subnetworks %D cannot be less than 1",Nsubnet); 58e2aaf10cSShri Abhyankar 59e2aaf10cSShri Abhyankar if(Nsubnet > 0) PetscValidLogicalCollectiveInt(dm,Nsubnet,2); 60e2aaf10cSShri Abhyankar if(network->nsubnet != 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Network sizes alread set, cannot resize the network"); 61e2aaf10cSShri Abhyankar 62e2aaf10cSShri Abhyankar network->nsubnet = Nsubnet; 63e2aaf10cSShri Abhyankar ierr = PetscCalloc1(Nsubnet,&network->subnet);CHKERRQ(ierr); 64e2aaf10cSShri Abhyankar for(i=0; i < network->nsubnet; i++) { 65e2aaf10cSShri Abhyankar if (NV[i] > 0) PetscValidLogicalCollectiveInt(dm,NV[i],5); 66e2aaf10cSShri Abhyankar if (NE[i] > 0) PetscValidLogicalCollectiveInt(dm,NE[i],6); 67e2aaf10cSShri Abhyankar if (NV[i] > 0 && nV[i] > NV[i]) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Subnetwork %D: Local vertex size %D cannot be larger than global vertex size %D",i,nV[i],NV[i]); 68e2aaf10cSShri Abhyankar if (NE[i] > 0 && nE[i] > NE[i]) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Subnetwork %D: Local edge size %D cannot be larger than global edge size %D",i,nE[i],NE[i]); 69e2aaf10cSShri Abhyankar a[0] = nV[i]; a[1] = nE[i]; 70b2566f29SBarry Smith ierr = MPIU_Allreduce(a,b,2,MPIU_INT,MPI_SUM,PetscObjectComm((PetscObject)dm));CHKERRQ(ierr); 71e2aaf10cSShri Abhyankar network->subnet[i].Nvtx = b[0]; network->subnet[i].Nedge = b[1]; 72e2aaf10cSShri Abhyankar 73e2aaf10cSShri Abhyankar network->subnet[i].id = i; 74e2aaf10cSShri Abhyankar 75e2aaf10cSShri Abhyankar network->subnet[i].nvtx = nV[i]; 76e2aaf10cSShri Abhyankar network->subnet[i].vStart = network->nVertices; 77e2aaf10cSShri Abhyankar network->subnet[i].vEnd = network->subnet[i].vStart + network->subnet[i].Nvtx; 78e2aaf10cSShri Abhyankar network->nVertices += network->subnet[i].nvtx; 79e2aaf10cSShri Abhyankar network->NVertices += network->subnet[i].Nvtx; 80e2aaf10cSShri Abhyankar 81e2aaf10cSShri Abhyankar network->subnet[i].nedge = nE[i]; 82e2aaf10cSShri Abhyankar network->subnet[i].eStart = network->nEdges; 83e2aaf10cSShri Abhyankar network->subnet[i].eEnd = network->subnet[i].eStart + network->subnet[i].Nedge; 84e2aaf10cSShri Abhyankar network->nEdges += network->subnet[i].nedge; 85e2aaf10cSShri Abhyankar network->NEdges += network->subnet[i].Nedge; 865f2c45f1SShri Abhyankar } 875f2c45f1SShri Abhyankar PetscFunctionReturn(0); 885f2c45f1SShri Abhyankar } 895f2c45f1SShri Abhyankar 905f2c45f1SShri Abhyankar /*@ 915f2c45f1SShri Abhyankar DMNetworkSetEdgeList - Sets the list of local edges (vertex connectivity) for the network 925f2c45f1SShri Abhyankar 935f2c45f1SShri Abhyankar Logically collective on DM 945f2c45f1SShri Abhyankar 955f2c45f1SShri Abhyankar Input Parameters: 96e2aaf10cSShri Abhyankar . edges - list of edges for each subnetwork 975f2c45f1SShri Abhyankar 985f2c45f1SShri Abhyankar Notes: 995f2c45f1SShri Abhyankar There is no copy involved in this operation, only the pointer is referenced. The edgelist should 1005f2c45f1SShri Abhyankar not be destroyed before the call to DMNetworkLayoutSetUp 1015f2c45f1SShri Abhyankar 1025f2c45f1SShri Abhyankar Level: intermediate 1035f2c45f1SShri Abhyankar 1045f2c45f1SShri Abhyankar .seealso: DMNetworkCreate, DMNetworkSetSizes 1055f2c45f1SShri Abhyankar @*/ 106e2aaf10cSShri Abhyankar PetscErrorCode DMNetworkSetEdgeList(DM dm, int *edgelist[]) 1075f2c45f1SShri Abhyankar { 1085f2c45f1SShri Abhyankar DM_Network *network = (DM_Network*) dm->data; 109e2aaf10cSShri Abhyankar PetscInt i; 1105f2c45f1SShri Abhyankar 1115f2c45f1SShri Abhyankar PetscFunctionBegin; 112e2aaf10cSShri Abhyankar for(i=0; i < network->nsubnet; i++) { 113e2aaf10cSShri Abhyankar network->subnet[i].edgelist = edgelist[i]; 114e2aaf10cSShri Abhyankar } 1155f2c45f1SShri Abhyankar PetscFunctionReturn(0); 1165f2c45f1SShri Abhyankar } 1175f2c45f1SShri Abhyankar 1185f2c45f1SShri Abhyankar /*@ 1195f2c45f1SShri Abhyankar DMNetworkLayoutSetUp - Sets up the bare layout (graph) for the network 1205f2c45f1SShri Abhyankar 1215f2c45f1SShri Abhyankar Collective on DM 1225f2c45f1SShri Abhyankar 1235f2c45f1SShri Abhyankar Input Parameters 1245f2c45f1SShri Abhyankar . DM - the dmnetwork object 1255f2c45f1SShri Abhyankar 1265f2c45f1SShri Abhyankar Notes: 1275f2c45f1SShri Abhyankar This routine should be called after the network sizes and edgelists have been provided. It creates 1285f2c45f1SShri Abhyankar the bare layout of the network and sets up the network to begin insertion of components. 1295f2c45f1SShri Abhyankar 1305f2c45f1SShri Abhyankar All the components should be registered before calling this routine. 1315f2c45f1SShri Abhyankar 1325f2c45f1SShri Abhyankar Level: intermediate 1335f2c45f1SShri Abhyankar 1345f2c45f1SShri Abhyankar .seealso: DMNetworkSetSizes, DMNetworkSetEdgeList 1355f2c45f1SShri Abhyankar @*/ 1365f2c45f1SShri Abhyankar PetscErrorCode DMNetworkLayoutSetUp(DM dm) 1375f2c45f1SShri Abhyankar { 1385f2c45f1SShri Abhyankar PetscErrorCode ierr; 1395f2c45f1SShri Abhyankar DM_Network *network = (DM_Network*) dm->data; 1405f2c45f1SShri Abhyankar PetscInt dim = 1; /* One dimensional network */ 1415f2c45f1SShri Abhyankar PetscInt numCorners=2; 1425f2c45f1SShri Abhyankar PetscInt spacedim=2; 1435f2c45f1SShri Abhyankar double *vertexcoords=NULL; 144e2aaf10cSShri Abhyankar PetscInt i,j; 1455f2c45f1SShri Abhyankar PetscInt ndata; 146e2aaf10cSShri Abhyankar PetscInt ctr=0; 1475f2c45f1SShri Abhyankar 1485f2c45f1SShri Abhyankar PetscFunctionBegin; 149e2aaf10cSShri Abhyankar 1506fefedf4SHong Zhang if (network->nVertices) { 1516fefedf4SHong Zhang ierr = PetscCalloc1(numCorners*network->nVertices,&vertexcoords);CHKERRQ(ierr); 1525f2c45f1SShri Abhyankar } 153e2aaf10cSShri Abhyankar 154e2aaf10cSShri Abhyankar /* Create the edgelist for the network by concatenating edgelists of the subnetworks */ 155e2aaf10cSShri Abhyankar ierr = PetscCalloc1(2*network->nEdges,&network->edges);CHKERRQ(ierr); 156e2aaf10cSShri Abhyankar for(i=0; i < network->nsubnet; i++) { 157e2aaf10cSShri Abhyankar for(j = 0; j < network->subnet[i].nedge; j++) { 158e2aaf10cSShri Abhyankar network->edges[2*ctr] = network->subnet[i].vStart + network->subnet[i].edgelist[2*j]; 159e2aaf10cSShri Abhyankar network->edges[2*ctr+1] = network->subnet[i].vStart + network->subnet[i].edgelist[2*j+1]; 160e2aaf10cSShri Abhyankar ctr++; 161e2aaf10cSShri Abhyankar } 162e2aaf10cSShri Abhyankar } 163e2aaf10cSShri Abhyankar 164e2aaf10cSShri Abhyankar #if 0 165e2aaf10cSShri Abhyankar for(i=0; i < network->nEdges; i++) { 166e2aaf10cSShri Abhyankar ierr = PetscPrintf(PETSC_COMM_SELF,"[%D %D]",network->edges[2*i],network->edges[2*i+1]);CHKERRQ(ierr); 167e2aaf10cSShri Abhyankar } 168e2aaf10cSShri Abhyankar #endif 169e2aaf10cSShri Abhyankar 1706fefedf4SHong Zhang ierr = DMPlexCreateFromCellList(PetscObjectComm((PetscObject)dm),dim,network->nEdges,network->nVertices,numCorners,PETSC_FALSE,network->edges,spacedim,vertexcoords,&network->plex);CHKERRQ(ierr); 1716fefedf4SHong Zhang if (network->nVertices) { 1725f2c45f1SShri Abhyankar ierr = PetscFree(vertexcoords);CHKERRQ(ierr); 1735f2c45f1SShri Abhyankar } 174e2aaf10cSShri Abhyankar ierr = PetscFree(network->edges);CHKERRQ(ierr); 175e2aaf10cSShri Abhyankar 1765f2c45f1SShri Abhyankar ierr = DMPlexGetChart(network->plex,&network->pStart,&network->pEnd);CHKERRQ(ierr); 1775f2c45f1SShri Abhyankar ierr = DMPlexGetHeightStratum(network->plex,0,&network->eStart,&network->eEnd);CHKERRQ(ierr); 1785f2c45f1SShri Abhyankar ierr = DMPlexGetHeightStratum(network->plex,1,&network->vStart,&network->vEnd);CHKERRQ(ierr); 1795f2c45f1SShri Abhyankar 1805f2c45f1SShri Abhyankar ierr = PetscSectionCreate(PetscObjectComm((PetscObject)dm),&network->DataSection);CHKERRQ(ierr); 1815f2c45f1SShri Abhyankar ierr = PetscSectionCreate(PetscObjectComm((PetscObject)dm),&network->DofSection);CHKERRQ(ierr); 1825f2c45f1SShri Abhyankar ierr = PetscSectionSetChart(network->DataSection,network->pStart,network->pEnd);CHKERRQ(ierr); 1835f2c45f1SShri Abhyankar ierr = PetscSectionSetChart(network->DofSection,network->pStart,network->pEnd);CHKERRQ(ierr); 1845f2c45f1SShri Abhyankar 185*2727e31bSShri Abhyankar /* Create vertices and edges array for the subnetworks */ 186*2727e31bSShri Abhyankar for(j=0; j < network->nsubnet; j++) { 187*2727e31bSShri Abhyankar ierr = PetscCalloc1(network->subnet[j].nedge,&network->subnet[j].edges);CHKERRQ(ierr); 188*2727e31bSShri Abhyankar ierr = PetscCalloc1(network->subnet[j].nvtx,&network->subnet[j].vertices);CHKERRQ(ierr); 189*2727e31bSShri Abhyankar /* Temporarily setting nvtx and nedge to 0 so we can use them as counters in the below for loop. 190*2727e31bSShri Abhyankar These get updated when the vertices and edges are added. */ 191*2727e31bSShri Abhyankar network->subnet[j].nvtx = network->subnet[j].nedge = 0; 192*2727e31bSShri Abhyankar 193*2727e31bSShri Abhyankar } 194*2727e31bSShri Abhyankar 1955f2c45f1SShri Abhyankar network->dataheadersize = sizeof(struct _p_DMNetworkComponentHeader)/sizeof(DMNetworkComponentGenericDataType); 1966caa05f4SBarry Smith ierr = PetscCalloc1(network->pEnd-network->pStart,&network->header);CHKERRQ(ierr); 197e2aaf10cSShri Abhyankar for(i=network->eStart; i < network->eEnd; i++) { 198e2aaf10cSShri Abhyankar network->header[i].index = i; /* Global edge number */ 199e2aaf10cSShri Abhyankar for(j=0; j < network->nsubnet; j++) { 200e2aaf10cSShri Abhyankar if((network->subnet[j].eStart <= i) && (i < network->subnet[j].eEnd)) { 201e2aaf10cSShri Abhyankar network->header[i].subnetid = j; /* Subnetwork id */ 202*2727e31bSShri Abhyankar network->subnet[j].edges[network->subnet[j].nedge++] = i; 203e2aaf10cSShri Abhyankar break; 204e2aaf10cSShri Abhyankar } 2057b6afd5bSHong Zhang } 2065f2c45f1SShri Abhyankar network->header[i].ndata = 0; 2075f2c45f1SShri Abhyankar ndata = network->header[i].ndata; 2085f2c45f1SShri Abhyankar ierr = PetscSectionAddDof(network->DataSection,i,network->dataheadersize);CHKERRQ(ierr); 2095f2c45f1SShri Abhyankar network->header[i].offset[ndata] = 0; 2105f2c45f1SShri Abhyankar } 211e2aaf10cSShri Abhyankar 212e2aaf10cSShri Abhyankar for(i=network->vStart; i < network->vEnd; i++) { 213e2aaf10cSShri Abhyankar network->header[i].index = i - network->vStart; 214e2aaf10cSShri Abhyankar for(j=0; j < network->nsubnet; j++) { 215e2aaf10cSShri Abhyankar if((network->subnet[j].vStart <= i-network->vStart) && (i-network->vStart < network->subnet[j].vEnd)) { 216e2aaf10cSShri Abhyankar network->header[i].subnetid = j; 217*2727e31bSShri Abhyankar network->subnet[j].vertices[network->subnet[j].nvtx++] = i; 218e2aaf10cSShri Abhyankar break; 219e2aaf10cSShri Abhyankar } 220e2aaf10cSShri Abhyankar } 221e2aaf10cSShri Abhyankar network->header[i].ndata = 0; 222e2aaf10cSShri Abhyankar ndata = network->header[i].ndata; 223e2aaf10cSShri Abhyankar ierr = PetscSectionAddDof(network->DataSection,i,network->dataheadersize);CHKERRQ(ierr); 224e2aaf10cSShri Abhyankar network->header[i].offset[ndata] = 0; 225e2aaf10cSShri Abhyankar } 226e2aaf10cSShri Abhyankar 227854ce69bSBarry Smith ierr = PetscMalloc1(network->pEnd-network->pStart,&network->cvalue);CHKERRQ(ierr); 228e2aaf10cSShri Abhyankar 2295f2c45f1SShri Abhyankar PetscFunctionReturn(0); 2305f2c45f1SShri Abhyankar } 2315f2c45f1SShri Abhyankar 23294ef8ddeSSatish Balay /*@C 233*2727e31bSShri Abhyankar DMNetworkGetSubnetworkInfo - Returns the info for the subnetwork 234*2727e31bSShri Abhyankar 235*2727e31bSShri Abhyankar Input Parameters 236*2727e31bSShri Abhyankar + dm - the number object 237*2727e31bSShri Abhyankar - id - the ID (integer) of the subnetwork 238*2727e31bSShri Abhyankar 239*2727e31bSShri Abhyankar Output Parameters 240*2727e31bSShri Abhyankar + nv - number of vertices (local) 241*2727e31bSShri Abhyankar . ne - number of edges (local) 242*2727e31bSShri Abhyankar . vtx - local vertices for this subnetwork 243*2727e31bSShri Abhyankar . edge - local edges for this subnetwork 244*2727e31bSShri Abhyankar 245*2727e31bSShri Abhyankar Notes: 246*2727e31bSShri Abhyankar Cannot call this routine before DMNetworkLayoutSetup() 247*2727e31bSShri Abhyankar 248*2727e31bSShri Abhyankar .seealso: DMNetworkLayoutSetUp, DMNetworkCreate 249*2727e31bSShri Abhyankar @*/ 250*2727e31bSShri Abhyankar PetscErrorCode DMNetworkGetSubnetworkInfo(DM netdm,PetscInt id,PetscInt *nv, PetscInt *ne,const PetscInt **vtx, const PetscInt **edge) 251*2727e31bSShri Abhyankar { 252*2727e31bSShri Abhyankar DM_Network *network = (DM_Network*) netdm->data; 253*2727e31bSShri Abhyankar 254*2727e31bSShri Abhyankar PetscFunctionBegin; 255*2727e31bSShri Abhyankar *nv = network->subnet[id].nvtx; 256*2727e31bSShri Abhyankar *ne = network->subnet[id].nedge; 257*2727e31bSShri Abhyankar *vtx = network->subnet[id].vertices; 258*2727e31bSShri Abhyankar *edge = network->subnet[id].edges; 259*2727e31bSShri Abhyankar PetscFunctionReturn(0); 260*2727e31bSShri Abhyankar } 261*2727e31bSShri Abhyankar 262*2727e31bSShri Abhyankar /*@C 2635f2c45f1SShri Abhyankar DMNetworkRegisterComponent - Registers the network component 2645f2c45f1SShri Abhyankar 2655f2c45f1SShri Abhyankar Logically collective on DM 2665f2c45f1SShri Abhyankar 2675f2c45f1SShri Abhyankar Input Parameters 2685f2c45f1SShri Abhyankar + dm - the network object 2695f2c45f1SShri Abhyankar . name - the component name 2705f2c45f1SShri Abhyankar - size - the storage size in bytes for this component data 2715f2c45f1SShri Abhyankar 2725f2c45f1SShri Abhyankar Output Parameters 2735f2c45f1SShri Abhyankar . key - an integer key that defines the component 2745f2c45f1SShri Abhyankar 2755f2c45f1SShri Abhyankar Notes 2765f2c45f1SShri Abhyankar This routine should be called by all processors before calling DMNetworkLayoutSetup(). 2775f2c45f1SShri Abhyankar 2785f2c45f1SShri Abhyankar Level: intermediate 2795f2c45f1SShri Abhyankar 2805f2c45f1SShri Abhyankar .seealso: DMNetworkLayoutSetUp, DMNetworkCreate 2815f2c45f1SShri Abhyankar @*/ 2825f2c45f1SShri Abhyankar PetscErrorCode DMNetworkRegisterComponent(DM dm,const char *name,PetscInt size,PetscInt *key) 2835f2c45f1SShri Abhyankar { 2845f2c45f1SShri Abhyankar PetscErrorCode ierr; 2855f2c45f1SShri Abhyankar DM_Network *network = (DM_Network*) dm->data; 2865f2c45f1SShri Abhyankar DMNetworkComponent *component=&network->component[network->ncomponent]; 2875f2c45f1SShri Abhyankar PetscBool flg=PETSC_FALSE; 2885f2c45f1SShri Abhyankar PetscInt i; 2895f2c45f1SShri Abhyankar 2905f2c45f1SShri Abhyankar PetscFunctionBegin; 2915f2c45f1SShri Abhyankar 2925f2c45f1SShri Abhyankar for (i=0; i < network->ncomponent; i++) { 2935f2c45f1SShri Abhyankar ierr = PetscStrcmp(component->name,name,&flg);CHKERRQ(ierr); 2945f2c45f1SShri Abhyankar if (flg) { 2955f2c45f1SShri Abhyankar *key = i; 2965f2c45f1SShri Abhyankar PetscFunctionReturn(0); 2975f2c45f1SShri Abhyankar } 2985f2c45f1SShri Abhyankar } 2995f2c45f1SShri Abhyankar 3005f2c45f1SShri Abhyankar ierr = PetscStrcpy(component->name,name);CHKERRQ(ierr); 3015f2c45f1SShri Abhyankar component->size = size/sizeof(DMNetworkComponentGenericDataType); 3025f2c45f1SShri Abhyankar *key = network->ncomponent; 3035f2c45f1SShri Abhyankar network->ncomponent++; 3045f2c45f1SShri Abhyankar PetscFunctionReturn(0); 3055f2c45f1SShri Abhyankar } 3065f2c45f1SShri Abhyankar 3075f2c45f1SShri Abhyankar /*@ 3085f2c45f1SShri Abhyankar DMNetworkGetVertexRange - Get the bounds [start, end) for the vertices. 3095f2c45f1SShri Abhyankar 3105f2c45f1SShri Abhyankar Not Collective 3115f2c45f1SShri Abhyankar 3125f2c45f1SShri Abhyankar Input Parameters: 3135f2c45f1SShri Abhyankar + dm - The DMNetwork object 3145f2c45f1SShri Abhyankar 3155f2c45f1SShri Abhyankar Output Paramters: 3165f2c45f1SShri Abhyankar + vStart - The first vertex point 3175f2c45f1SShri Abhyankar - vEnd - One beyond the last vertex point 3185f2c45f1SShri Abhyankar 3195f2c45f1SShri Abhyankar Level: intermediate 3205f2c45f1SShri Abhyankar 3215f2c45f1SShri Abhyankar .seealso: DMNetworkGetEdgeRange 3225f2c45f1SShri Abhyankar @*/ 3235f2c45f1SShri Abhyankar PetscErrorCode DMNetworkGetVertexRange(DM dm,PetscInt *vStart,PetscInt *vEnd) 3245f2c45f1SShri Abhyankar { 3255f2c45f1SShri Abhyankar DM_Network *network = (DM_Network*)dm->data; 3265f2c45f1SShri Abhyankar 3275f2c45f1SShri Abhyankar PetscFunctionBegin; 3285f2c45f1SShri Abhyankar if (vStart) *vStart = network->vStart; 3295f2c45f1SShri Abhyankar if (vEnd) *vEnd = network->vEnd; 3305f2c45f1SShri Abhyankar PetscFunctionReturn(0); 3315f2c45f1SShri Abhyankar } 3325f2c45f1SShri Abhyankar 3335f2c45f1SShri Abhyankar /*@ 3345f2c45f1SShri Abhyankar DMNetworkGetEdgeRange - Get the bounds [start, end) for the edges. 3355f2c45f1SShri Abhyankar 3365f2c45f1SShri Abhyankar Not Collective 3375f2c45f1SShri Abhyankar 3385f2c45f1SShri Abhyankar Input Parameters: 3395f2c45f1SShri Abhyankar + dm - The DMNetwork object 3405f2c45f1SShri Abhyankar 3415f2c45f1SShri Abhyankar Output Paramters: 3425f2c45f1SShri Abhyankar + eStart - The first edge point 3435f2c45f1SShri Abhyankar - eEnd - One beyond the last edge point 3445f2c45f1SShri Abhyankar 3455f2c45f1SShri Abhyankar Level: intermediate 3465f2c45f1SShri Abhyankar 3475f2c45f1SShri Abhyankar .seealso: DMNetworkGetVertexRange 3485f2c45f1SShri Abhyankar @*/ 3495f2c45f1SShri Abhyankar PetscErrorCode DMNetworkGetEdgeRange(DM dm,PetscInt *eStart,PetscInt *eEnd) 3505f2c45f1SShri Abhyankar { 3515f2c45f1SShri Abhyankar DM_Network *network = (DM_Network*)dm->data; 3525f2c45f1SShri Abhyankar 3535f2c45f1SShri Abhyankar PetscFunctionBegin; 3545f2c45f1SShri Abhyankar if (eStart) *eStart = network->eStart; 3555f2c45f1SShri Abhyankar if (eEnd) *eEnd = network->eEnd; 3565f2c45f1SShri Abhyankar PetscFunctionReturn(0); 3575f2c45f1SShri Abhyankar } 3585f2c45f1SShri Abhyankar 3597b6afd5bSHong Zhang /*@ 360e85e6aecSHong Zhang DMNetworkGetGlobalEdgeIndex - Get the user global numbering for the edge. 3617b6afd5bSHong Zhang 3627b6afd5bSHong Zhang Not Collective 3637b6afd5bSHong Zhang 3647b6afd5bSHong Zhang Input Parameters: 3657b6afd5bSHong Zhang + dm - DMNetwork object 366e85e6aecSHong Zhang - p - edge point 3677b6afd5bSHong Zhang 3687b6afd5bSHong Zhang Output Paramters: 369e85e6aecSHong Zhang . index - user global numbering for the edge 3707b6afd5bSHong Zhang 3717b6afd5bSHong Zhang Level: intermediate 3727b6afd5bSHong Zhang 373e85e6aecSHong Zhang .seealso: DMNetworkGetGlobalVertexIndex 3747b6afd5bSHong Zhang @*/ 375e85e6aecSHong Zhang PetscErrorCode DMNetworkGetGlobalEdgeIndex(DM dm,PetscInt p,PetscInt *index) 3767b6afd5bSHong Zhang { 3777b6afd5bSHong Zhang PetscErrorCode ierr; 3787b6afd5bSHong Zhang DM_Network *network = (DM_Network*)dm->data; 3797b6afd5bSHong Zhang PetscInt offsetp; 3807b6afd5bSHong Zhang DMNetworkComponentHeader header; 3817b6afd5bSHong Zhang 3827b6afd5bSHong Zhang PetscFunctionBegin; 3837b6afd5bSHong Zhang ierr = PetscSectionGetOffset(network->DataSection,p,&offsetp);CHKERRQ(ierr); 3847b6afd5bSHong Zhang header = (DMNetworkComponentHeader)(network->componentdataarray+offsetp); 385e85e6aecSHong Zhang *index = header->index; 3867b6afd5bSHong Zhang PetscFunctionReturn(0); 3877b6afd5bSHong Zhang } 3887b6afd5bSHong Zhang 3895f2c45f1SShri Abhyankar /*@ 390e85e6aecSHong Zhang DMNetworkGetGlobalVertexIndex - Get the user global numbering for the vertex. 391e85e6aecSHong Zhang 392e85e6aecSHong Zhang Not Collective 393e85e6aecSHong Zhang 394e85e6aecSHong Zhang Input Parameters: 395e85e6aecSHong Zhang + dm - DMNetwork object 396e85e6aecSHong Zhang - p - vertex point 397e85e6aecSHong Zhang 398e85e6aecSHong Zhang Output Paramters: 399e85e6aecSHong Zhang . index - user global numbering for the vertex 400e85e6aecSHong Zhang 401e85e6aecSHong Zhang Level: intermediate 402e85e6aecSHong Zhang 403e85e6aecSHong Zhang .seealso: DMNetworkGetGlobalEdgeIndex 404e85e6aecSHong Zhang @*/ 405e85e6aecSHong Zhang PetscErrorCode DMNetworkGetGlobalVertexIndex(DM dm,PetscInt p,PetscInt *index) 406e85e6aecSHong Zhang { 407e85e6aecSHong Zhang PetscErrorCode ierr; 408e85e6aecSHong Zhang DM_Network *network = (DM_Network*)dm->data; 409e85e6aecSHong Zhang PetscInt offsetp; 410e85e6aecSHong Zhang DMNetworkComponentHeader header; 411e85e6aecSHong Zhang 412e85e6aecSHong Zhang PetscFunctionBegin; 413e85e6aecSHong Zhang ierr = PetscSectionGetOffset(network->DataSection,p,&offsetp);CHKERRQ(ierr); 414e85e6aecSHong Zhang header = (DMNetworkComponentHeader)(network->componentdataarray+offsetp); 415e85e6aecSHong Zhang *index = header->index; 416e85e6aecSHong Zhang PetscFunctionReturn(0); 417e85e6aecSHong Zhang } 418e85e6aecSHong Zhang 419e85e6aecSHong Zhang /*@ 420325661f6SSatish Balay DMNetworkAddComponent - Adds a network component at the given point (vertex/edge) 4215f2c45f1SShri Abhyankar 4225f2c45f1SShri Abhyankar Not Collective 4235f2c45f1SShri Abhyankar 4245f2c45f1SShri Abhyankar Input Parameters: 4255f2c45f1SShri Abhyankar + dm - The DMNetwork object 4265f2c45f1SShri Abhyankar . p - vertex/edge point 4275f2c45f1SShri Abhyankar . componentkey - component key returned while registering the component 4285f2c45f1SShri Abhyankar - compvalue - pointer to the data structure for the component 4295f2c45f1SShri Abhyankar 4305f2c45f1SShri Abhyankar Level: intermediate 4315f2c45f1SShri Abhyankar 4325f2c45f1SShri Abhyankar .seealso: DMNetworkGetVertexRange, DMNetworkGetEdgeRange, DMNetworkRegisterComponent 4335f2c45f1SShri Abhyankar @*/ 4345f2c45f1SShri Abhyankar PetscErrorCode DMNetworkAddComponent(DM dm, PetscInt p,PetscInt componentkey,void* compvalue) 4355f2c45f1SShri Abhyankar { 4365f2c45f1SShri Abhyankar DM_Network *network = (DM_Network*)dm->data; 43743a39a44SBarry Smith DMNetworkComponent *component = &network->component[componentkey]; 4385f2c45f1SShri Abhyankar DMNetworkComponentHeader header = &network->header[p]; 4395f2c45f1SShri Abhyankar DMNetworkComponentValue cvalue = &network->cvalue[p]; 4405f2c45f1SShri Abhyankar PetscErrorCode ierr; 4415f2c45f1SShri Abhyankar 4425f2c45f1SShri Abhyankar PetscFunctionBegin; 443fa58f0a9SHong Zhang if (header->ndata == MAX_DATA_AT_POINT) SETERRQ1(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_OUTOFRANGE,"Number of components at a point exceeds the max %D",MAX_DATA_AT_POINT); 444fa58f0a9SHong Zhang 44543a39a44SBarry Smith header->size[header->ndata] = component->size; 44643a39a44SBarry Smith ierr = PetscSectionAddDof(network->DataSection,p,component->size);CHKERRQ(ierr); 4475f2c45f1SShri Abhyankar header->key[header->ndata] = componentkey; 4485f2c45f1SShri Abhyankar if (header->ndata != 0) header->offset[header->ndata] = header->offset[header->ndata-1] + header->size[header->ndata-1]; 4495f2c45f1SShri Abhyankar 4505f2c45f1SShri Abhyankar cvalue->data[header->ndata] = (void*)compvalue; 4515f2c45f1SShri Abhyankar header->ndata++; 4525f2c45f1SShri Abhyankar PetscFunctionReturn(0); 4535f2c45f1SShri Abhyankar } 4545f2c45f1SShri Abhyankar 4555f2c45f1SShri Abhyankar /*@ 4565f2c45f1SShri Abhyankar DMNetworkGetNumComponents - Get the number of components at a vertex/edge 4575f2c45f1SShri Abhyankar 4585f2c45f1SShri Abhyankar Not Collective 4595f2c45f1SShri Abhyankar 4605f2c45f1SShri Abhyankar Input Parameters: 4615f2c45f1SShri Abhyankar + dm - The DMNetwork object 4625f2c45f1SShri Abhyankar . p - vertex/edge point 4635f2c45f1SShri Abhyankar 4645f2c45f1SShri Abhyankar Output Parameters: 4655f2c45f1SShri Abhyankar . numcomponents - Number of components at the vertex/edge 4665f2c45f1SShri Abhyankar 4675f2c45f1SShri Abhyankar Level: intermediate 4685f2c45f1SShri Abhyankar 4695f2c45f1SShri Abhyankar .seealso: DMNetworkRegisterComponent, DMNetworkAddComponent 4705f2c45f1SShri Abhyankar @*/ 4715f2c45f1SShri Abhyankar PetscErrorCode DMNetworkGetNumComponents(DM dm,PetscInt p,PetscInt *numcomponents) 4725f2c45f1SShri Abhyankar { 4735f2c45f1SShri Abhyankar PetscErrorCode ierr; 4745f2c45f1SShri Abhyankar PetscInt offset; 4755f2c45f1SShri Abhyankar DM_Network *network = (DM_Network*)dm->data; 4765f2c45f1SShri Abhyankar 4775f2c45f1SShri Abhyankar PetscFunctionBegin; 4785f2c45f1SShri Abhyankar ierr = PetscSectionGetOffset(network->DataSection,p,&offset);CHKERRQ(ierr); 4795f2c45f1SShri Abhyankar *numcomponents = ((DMNetworkComponentHeader)(network->componentdataarray+offset))->ndata; 4805f2c45f1SShri Abhyankar PetscFunctionReturn(0); 4815f2c45f1SShri Abhyankar } 4825f2c45f1SShri Abhyankar 4835f2c45f1SShri Abhyankar /*@ 484a730d845SHong Zhang DMNetworkGetComponentKeyOffset - Gets the type along with the offset for indexing the 4855f2c45f1SShri Abhyankar component value from the component data array 4865f2c45f1SShri Abhyankar 4875f2c45f1SShri Abhyankar Not Collective 4885f2c45f1SShri Abhyankar 4895f2c45f1SShri Abhyankar Input Parameters: 4905f2c45f1SShri Abhyankar + dm - The DMNetwork object 4915f2c45f1SShri Abhyankar . p - vertex/edge point 4925f2c45f1SShri Abhyankar - compnum - component number 4935f2c45f1SShri Abhyankar 4945f2c45f1SShri Abhyankar Output Parameters: 4955f2c45f1SShri Abhyankar + compkey - the key obtained when registering the component 4965f2c45f1SShri Abhyankar - offset - offset into the component data array associated with the vertex/edge point 4975f2c45f1SShri Abhyankar 4985f2c45f1SShri Abhyankar Notes: 4995f2c45f1SShri Abhyankar Typical usage: 5005f2c45f1SShri Abhyankar 5015f2c45f1SShri Abhyankar DMNetworkGetComponentDataArray(dm, &arr); 5025f2c45f1SShri Abhyankar DMNetworkGetVertex/EdgeRange(dm,&Start,&End); 5035f2c45f1SShri Abhyankar Loop over vertices or edges 5045f2c45f1SShri Abhyankar DMNetworkGetNumComponents(dm,v,&numcomps); 5055f2c45f1SShri Abhyankar Loop over numcomps 506a730d845SHong Zhang DMNetworkGetComponentKeyOffset(dm,v,compnum,&key,&offset); 5075f2c45f1SShri Abhyankar compdata = (UserCompDataType)(arr+offset); 5085f2c45f1SShri Abhyankar 5095f2c45f1SShri Abhyankar Level: intermediate 5105f2c45f1SShri Abhyankar 5115f2c45f1SShri Abhyankar .seealso: DMNetworkGetNumComponents, DMNetworkGetComponentDataArray, 5125f2c45f1SShri Abhyankar @*/ 513a730d845SHong Zhang PetscErrorCode DMNetworkGetComponentKeyOffset(DM dm,PetscInt p, PetscInt compnum, PetscInt *compkey, PetscInt *offset) 5145f2c45f1SShri Abhyankar { 5155f2c45f1SShri Abhyankar PetscErrorCode ierr; 5165f2c45f1SShri Abhyankar PetscInt offsetp; 5175f2c45f1SShri Abhyankar DMNetworkComponentHeader header; 5185f2c45f1SShri Abhyankar DM_Network *network = (DM_Network*)dm->data; 5195f2c45f1SShri Abhyankar 5205f2c45f1SShri Abhyankar PetscFunctionBegin; 5215f2c45f1SShri Abhyankar ierr = PetscSectionGetOffset(network->DataSection,p,&offsetp);CHKERRQ(ierr); 5225f2c45f1SShri Abhyankar header = (DMNetworkComponentHeader)(network->componentdataarray+offsetp); 523d36f4e81SHong Zhang if (compkey) *compkey = header->key[compnum]; 524d36f4e81SHong Zhang if (offset) *offset = offsetp+network->dataheadersize+header->offset[compnum]; 5255f2c45f1SShri Abhyankar PetscFunctionReturn(0); 5265f2c45f1SShri Abhyankar } 5275f2c45f1SShri Abhyankar 5285f2c45f1SShri Abhyankar /*@ 5295f2c45f1SShri Abhyankar DMNetworkGetVariableOffset - Get the offset for accessing the variable associated with the given vertex/edge from the local vector. 5305f2c45f1SShri Abhyankar 5315f2c45f1SShri Abhyankar Not Collective 5325f2c45f1SShri Abhyankar 5335f2c45f1SShri Abhyankar Input Parameters: 5345f2c45f1SShri Abhyankar + dm - The DMNetwork object 5355f2c45f1SShri Abhyankar - p - the edge/vertex point 5365f2c45f1SShri Abhyankar 5375f2c45f1SShri Abhyankar Output Parameters: 5385f2c45f1SShri Abhyankar . offset - the offset 5395f2c45f1SShri Abhyankar 5405f2c45f1SShri Abhyankar Level: intermediate 5415f2c45f1SShri Abhyankar 5425f2c45f1SShri Abhyankar .seealso: DMNetworkGetVariableGlobalOffset, DMGetLocalVector 5435f2c45f1SShri Abhyankar @*/ 5445f2c45f1SShri Abhyankar PetscErrorCode DMNetworkGetVariableOffset(DM dm,PetscInt p,PetscInt *offset) 5455f2c45f1SShri Abhyankar { 5465f2c45f1SShri Abhyankar PetscErrorCode ierr; 5475f2c45f1SShri Abhyankar DM_Network *network = (DM_Network*)dm->data; 5485f2c45f1SShri Abhyankar 5495f2c45f1SShri Abhyankar PetscFunctionBegin; 5505f78ed8bSShri Abhyankar ierr = PetscSectionGetOffset(network->plex->defaultSection,p,offset);CHKERRQ(ierr); 5515f2c45f1SShri Abhyankar PetscFunctionReturn(0); 5525f2c45f1SShri Abhyankar } 5535f2c45f1SShri Abhyankar 5545f2c45f1SShri Abhyankar /*@ 5555f2c45f1SShri Abhyankar DMNetworkGetVariableGlobalOffset - Get the global offset for the variable associated with the given vertex/edge from the global vector. 5565f2c45f1SShri Abhyankar 5575f2c45f1SShri Abhyankar Not Collective 5585f2c45f1SShri Abhyankar 5595f2c45f1SShri Abhyankar Input Parameters: 5605f2c45f1SShri Abhyankar + dm - The DMNetwork object 5615f2c45f1SShri Abhyankar - p - the edge/vertex point 5625f2c45f1SShri Abhyankar 5635f2c45f1SShri Abhyankar Output Parameters: 5645f2c45f1SShri Abhyankar . offsetg - the offset 5655f2c45f1SShri Abhyankar 5665f2c45f1SShri Abhyankar Level: intermediate 5675f2c45f1SShri Abhyankar 5685f2c45f1SShri Abhyankar .seealso: DMNetworkGetVariableOffset, DMGetLocalVector 5695f2c45f1SShri Abhyankar @*/ 5705f2c45f1SShri Abhyankar PetscErrorCode DMNetworkGetVariableGlobalOffset(DM dm,PetscInt p,PetscInt *offsetg) 5715f2c45f1SShri Abhyankar { 5725f2c45f1SShri Abhyankar PetscErrorCode ierr; 5735f2c45f1SShri Abhyankar DM_Network *network = (DM_Network*)dm->data; 5745f2c45f1SShri Abhyankar 5755f2c45f1SShri Abhyankar PetscFunctionBegin; 5765f78ed8bSShri Abhyankar ierr = PetscSectionGetOffset(network->plex->defaultGlobalSection,p,offsetg);CHKERRQ(ierr); 5776fefedf4SHong Zhang if (*offsetg < 0) *offsetg = -(*offsetg + 1); /* Convert to actual global offset for ghost vertex */ 5785f2c45f1SShri Abhyankar PetscFunctionReturn(0); 5795f2c45f1SShri Abhyankar } 5805f2c45f1SShri Abhyankar 58124121865SAdrian Maldonado /*@ 58224121865SAdrian Maldonado DMNetworkGetEdgeOffset - Get the offset for accessing the variable associated with the given edge from the local subvector. 58324121865SAdrian Maldonado 58424121865SAdrian Maldonado Not Collective 58524121865SAdrian Maldonado 58624121865SAdrian Maldonado Input Parameters: 58724121865SAdrian Maldonado + dm - The DMNetwork object 58824121865SAdrian Maldonado - p - the edge point 58924121865SAdrian Maldonado 59024121865SAdrian Maldonado Output Parameters: 59124121865SAdrian Maldonado . offset - the offset 59224121865SAdrian Maldonado 59324121865SAdrian Maldonado Level: intermediate 59424121865SAdrian Maldonado 59524121865SAdrian Maldonado .seealso: DMNetworkGetVariableGlobalOffset, DMGetLocalVector 59624121865SAdrian Maldonado @*/ 59724121865SAdrian Maldonado PetscErrorCode DMNetworkGetEdgeOffset(DM dm,PetscInt p,PetscInt *offset) 59824121865SAdrian Maldonado { 59924121865SAdrian Maldonado PetscErrorCode ierr; 60024121865SAdrian Maldonado DM_Network *network = (DM_Network*)dm->data; 60124121865SAdrian Maldonado 60224121865SAdrian Maldonado PetscFunctionBegin; 60324121865SAdrian Maldonado 60424121865SAdrian Maldonado ierr = PetscSectionGetOffset(network->edge.DofSection,p,offset);CHKERRQ(ierr); 60524121865SAdrian Maldonado PetscFunctionReturn(0); 60624121865SAdrian Maldonado } 60724121865SAdrian Maldonado 60824121865SAdrian Maldonado /*@ 60924121865SAdrian Maldonado DMNetworkGetVertexOffset - Get the offset for accessing the variable associated with the given vertex from the local subvector. 61024121865SAdrian Maldonado 61124121865SAdrian Maldonado Not Collective 61224121865SAdrian Maldonado 61324121865SAdrian Maldonado Input Parameters: 61424121865SAdrian Maldonado + dm - The DMNetwork object 61524121865SAdrian Maldonado - p - the vertex point 61624121865SAdrian Maldonado 61724121865SAdrian Maldonado Output Parameters: 61824121865SAdrian Maldonado . offset - the offset 61924121865SAdrian Maldonado 62024121865SAdrian Maldonado Level: intermediate 62124121865SAdrian Maldonado 62224121865SAdrian Maldonado .seealso: DMNetworkGetVariableGlobalOffset, DMGetLocalVector 62324121865SAdrian Maldonado @*/ 62424121865SAdrian Maldonado PetscErrorCode DMNetworkGetVertexOffset(DM dm,PetscInt p,PetscInt *offset) 62524121865SAdrian Maldonado { 62624121865SAdrian Maldonado PetscErrorCode ierr; 62724121865SAdrian Maldonado DM_Network *network = (DM_Network*)dm->data; 62824121865SAdrian Maldonado 62924121865SAdrian Maldonado PetscFunctionBegin; 63024121865SAdrian Maldonado 63124121865SAdrian Maldonado p -= network->vStart; 63224121865SAdrian Maldonado 63324121865SAdrian Maldonado ierr = PetscSectionGetOffset(network->vertex.DofSection,p,offset);CHKERRQ(ierr); 63424121865SAdrian Maldonado PetscFunctionReturn(0); 63524121865SAdrian Maldonado } 6365f2c45f1SShri Abhyankar /*@ 6375f2c45f1SShri Abhyankar DMNetworkAddNumVariables - Add number of variables associated with a given point. 6385f2c45f1SShri Abhyankar 6395f2c45f1SShri Abhyankar Not Collective 6405f2c45f1SShri Abhyankar 6415f2c45f1SShri Abhyankar Input Parameters: 6425f2c45f1SShri Abhyankar + dm - The DMNetworkObject 6435f2c45f1SShri Abhyankar . p - the vertex/edge point 6445f2c45f1SShri Abhyankar - nvar - number of additional variables 6455f2c45f1SShri Abhyankar 6465f2c45f1SShri Abhyankar Level: intermediate 6475f2c45f1SShri Abhyankar 6485f2c45f1SShri Abhyankar .seealso: DMNetworkSetNumVariables 6495f2c45f1SShri Abhyankar @*/ 6505f2c45f1SShri Abhyankar PetscErrorCode DMNetworkAddNumVariables(DM dm,PetscInt p,PetscInt nvar) 6515f2c45f1SShri Abhyankar { 6525f2c45f1SShri Abhyankar PetscErrorCode ierr; 6535f2c45f1SShri Abhyankar DM_Network *network = (DM_Network*)dm->data; 6545f2c45f1SShri Abhyankar 6555f2c45f1SShri Abhyankar PetscFunctionBegin; 6565f2c45f1SShri Abhyankar ierr = PetscSectionAddDof(network->DofSection,p,nvar);CHKERRQ(ierr); 6575f2c45f1SShri Abhyankar PetscFunctionReturn(0); 6585f2c45f1SShri Abhyankar } 6595f2c45f1SShri Abhyankar 66027f51fceSHong Zhang /*@ 66127f51fceSHong Zhang DMNetworkGetNumVariables - Gets number of variables for a vertex/edge point. 66227f51fceSHong Zhang 66327f51fceSHong Zhang Not Collective 66427f51fceSHong Zhang 66527f51fceSHong Zhang Input Parameters: 66627f51fceSHong Zhang + dm - The DMNetworkObject 66727f51fceSHong Zhang - p - the vertex/edge point 66827f51fceSHong Zhang 66927f51fceSHong Zhang Output Parameters: 67027f51fceSHong Zhang . nvar - number of variables 67127f51fceSHong Zhang 67227f51fceSHong Zhang Level: intermediate 67327f51fceSHong Zhang 67427f51fceSHong Zhang .seealso: DMNetworkAddNumVariables, DMNetworkSddNumVariables 67527f51fceSHong Zhang @*/ 67627f51fceSHong Zhang PetscErrorCode DMNetworkGetNumVariables(DM dm,PetscInt p,PetscInt *nvar) 67727f51fceSHong Zhang { 67827f51fceSHong Zhang PetscErrorCode ierr; 67927f51fceSHong Zhang DM_Network *network = (DM_Network*)dm->data; 68027f51fceSHong Zhang 68127f51fceSHong Zhang PetscFunctionBegin; 68227f51fceSHong Zhang ierr = PetscSectionGetDof(network->DofSection,p,nvar);CHKERRQ(ierr); 68327f51fceSHong Zhang PetscFunctionReturn(0); 68427f51fceSHong Zhang } 68527f51fceSHong Zhang 6865f2c45f1SShri Abhyankar /*@ 6875f2c45f1SShri Abhyankar DMNetworkSetNumVariables - Sets number of variables for a vertex/edge point. 6885f2c45f1SShri Abhyankar 6895f2c45f1SShri Abhyankar Not Collective 6905f2c45f1SShri Abhyankar 6915f2c45f1SShri Abhyankar Input Parameters: 6925f2c45f1SShri Abhyankar + dm - The DMNetworkObject 6935f2c45f1SShri Abhyankar . p - the vertex/edge point 6945f2c45f1SShri Abhyankar - nvar - number of variables 6955f2c45f1SShri Abhyankar 6965f2c45f1SShri Abhyankar Level: intermediate 6975f2c45f1SShri Abhyankar 6985f2c45f1SShri Abhyankar .seealso: DMNetworkAddNumVariables 6995f2c45f1SShri Abhyankar @*/ 7005f2c45f1SShri Abhyankar PetscErrorCode DMNetworkSetNumVariables(DM dm,PetscInt p,PetscInt nvar) 7015f2c45f1SShri Abhyankar { 7025f2c45f1SShri Abhyankar PetscErrorCode ierr; 7035f2c45f1SShri Abhyankar DM_Network *network = (DM_Network*)dm->data; 7045f2c45f1SShri Abhyankar 7055f2c45f1SShri Abhyankar PetscFunctionBegin; 7065f2c45f1SShri Abhyankar ierr = PetscSectionSetDof(network->DofSection,p,nvar);CHKERRQ(ierr); 7075f2c45f1SShri Abhyankar PetscFunctionReturn(0); 7085f2c45f1SShri Abhyankar } 7095f2c45f1SShri Abhyankar 7105f2c45f1SShri Abhyankar /* Sets up the array that holds the data for all components and its associated section. This 7115f2c45f1SShri Abhyankar function is called during DMSetUp() */ 7125f2c45f1SShri Abhyankar PetscErrorCode DMNetworkComponentSetUp(DM dm) 7135f2c45f1SShri Abhyankar { 7145f2c45f1SShri Abhyankar PetscErrorCode ierr; 7155f2c45f1SShri Abhyankar DM_Network *network = (DM_Network*)dm->data; 7165f2c45f1SShri Abhyankar PetscInt arr_size; 7175f2c45f1SShri Abhyankar PetscInt p,offset,offsetp; 7185f2c45f1SShri Abhyankar DMNetworkComponentHeader header; 7195f2c45f1SShri Abhyankar DMNetworkComponentValue cvalue; 7205f2c45f1SShri Abhyankar DMNetworkComponentGenericDataType *componentdataarray; 7215f2c45f1SShri Abhyankar PetscInt ncomp, i; 7225f2c45f1SShri Abhyankar 7235f2c45f1SShri Abhyankar PetscFunctionBegin; 7245f2c45f1SShri Abhyankar ierr = PetscSectionSetUp(network->DataSection);CHKERRQ(ierr); 7255f2c45f1SShri Abhyankar ierr = PetscSectionGetStorageSize(network->DataSection,&arr_size);CHKERRQ(ierr); 72675b160a0SShri Abhyankar ierr = PetscMalloc1(arr_size,&network->componentdataarray);CHKERRQ(ierr); 7275f2c45f1SShri Abhyankar componentdataarray = network->componentdataarray; 7285f2c45f1SShri Abhyankar for (p = network->pStart; p < network->pEnd; p++) { 7295f2c45f1SShri Abhyankar ierr = PetscSectionGetOffset(network->DataSection,p,&offsetp);CHKERRQ(ierr); 7305f2c45f1SShri Abhyankar /* Copy header */ 7315f2c45f1SShri Abhyankar header = &network->header[p]; 732302440fdSBarry Smith ierr = PetscMemcpy(componentdataarray+offsetp,header,network->dataheadersize*sizeof(DMNetworkComponentGenericDataType));CHKERRQ(ierr); 7335f2c45f1SShri Abhyankar /* Copy data */ 7345f2c45f1SShri Abhyankar cvalue = &network->cvalue[p]; 7355f2c45f1SShri Abhyankar ncomp = header->ndata; 7365f2c45f1SShri Abhyankar for (i = 0; i < ncomp; i++) { 7375f2c45f1SShri Abhyankar offset = offsetp + network->dataheadersize + header->offset[i]; 738302440fdSBarry Smith ierr = PetscMemcpy(componentdataarray+offset,cvalue->data[i],header->size[i]*sizeof(DMNetworkComponentGenericDataType));CHKERRQ(ierr); 7395f2c45f1SShri Abhyankar } 7405f2c45f1SShri Abhyankar } 7415f2c45f1SShri Abhyankar PetscFunctionReturn(0); 7425f2c45f1SShri Abhyankar } 7435f2c45f1SShri Abhyankar 7445f2c45f1SShri Abhyankar /* Sets up the section for dofs. This routine is called during DMSetUp() */ 7455f2c45f1SShri Abhyankar PetscErrorCode DMNetworkVariablesSetUp(DM dm) 7465f2c45f1SShri Abhyankar { 7475f2c45f1SShri Abhyankar PetscErrorCode ierr; 7485f2c45f1SShri Abhyankar DM_Network *network = (DM_Network*)dm->data; 7495f2c45f1SShri Abhyankar 7505f2c45f1SShri Abhyankar PetscFunctionBegin; 7515f2c45f1SShri Abhyankar ierr = PetscSectionSetUp(network->DofSection);CHKERRQ(ierr); 7525f2c45f1SShri Abhyankar PetscFunctionReturn(0); 7535f2c45f1SShri Abhyankar } 7545f2c45f1SShri Abhyankar 7555f2c45f1SShri Abhyankar /*@C 7565f2c45f1SShri Abhyankar DMNetworkGetComponentDataArray - Returns the component data array 7575f2c45f1SShri Abhyankar 7585f2c45f1SShri Abhyankar Not Collective 7595f2c45f1SShri Abhyankar 7605f2c45f1SShri Abhyankar Input Parameters: 7615f2c45f1SShri Abhyankar . dm - The DMNetwork Object 7625f2c45f1SShri Abhyankar 7635f2c45f1SShri Abhyankar Output Parameters: 7645f2c45f1SShri Abhyankar . componentdataarray - array that holds data for all components 7655f2c45f1SShri Abhyankar 7665f2c45f1SShri Abhyankar Level: intermediate 7675f2c45f1SShri Abhyankar 768a730d845SHong Zhang .seealso: DMNetworkGetComponentKeyOffset, DMNetworkGetNumComponents 7695f2c45f1SShri Abhyankar @*/ 7705f2c45f1SShri Abhyankar PetscErrorCode DMNetworkGetComponentDataArray(DM dm,DMNetworkComponentGenericDataType **componentdataarray) 7715f2c45f1SShri Abhyankar { 7725f2c45f1SShri Abhyankar DM_Network *network = (DM_Network*)dm->data; 7735f2c45f1SShri Abhyankar 7745f2c45f1SShri Abhyankar PetscFunctionBegin; 7755f2c45f1SShri Abhyankar *componentdataarray = network->componentdataarray; 7765f2c45f1SShri Abhyankar PetscFunctionReturn(0); 7775f2c45f1SShri Abhyankar } 7785f2c45f1SShri Abhyankar 77924121865SAdrian Maldonado /* Get a subsection from a range of points */ 78024121865SAdrian Maldonado PetscErrorCode DMNetworkGetSubSection_private(PetscSection master, PetscInt pstart, PetscInt pend,PetscSection *subsection) 78124121865SAdrian Maldonado { 78224121865SAdrian Maldonado PetscErrorCode ierr; 78324121865SAdrian Maldonado PetscInt i, nvar; 78424121865SAdrian Maldonado 78524121865SAdrian Maldonado PetscFunctionBegin; 78624121865SAdrian Maldonado ierr = PetscSectionCreate(PetscObjectComm((PetscObject)master), subsection);CHKERRQ(ierr); 78724121865SAdrian Maldonado ierr = PetscSectionSetChart(*subsection, 0, pend - pstart);CHKERRQ(ierr); 78824121865SAdrian Maldonado for (i = pstart; i < pend; i++) { 78924121865SAdrian Maldonado ierr = PetscSectionGetDof(master,i,&nvar);CHKERRQ(ierr); 79024121865SAdrian Maldonado ierr = PetscSectionSetDof(*subsection, i - pstart, nvar);CHKERRQ(ierr); 79124121865SAdrian Maldonado } 79224121865SAdrian Maldonado 79324121865SAdrian Maldonado ierr = PetscSectionSetUp(*subsection);CHKERRQ(ierr); 79424121865SAdrian Maldonado PetscFunctionReturn(0); 79524121865SAdrian Maldonado } 79624121865SAdrian Maldonado 79724121865SAdrian Maldonado /* Create a submap of points with a GlobalToLocal structure */ 79824121865SAdrian Maldonado PetscErrorCode DMNetworkSetSubMap_private(PetscInt pstart, PetscInt pend, ISLocalToGlobalMapping *map) 79924121865SAdrian Maldonado { 80024121865SAdrian Maldonado PetscErrorCode ierr; 80124121865SAdrian Maldonado PetscInt i, *subpoints; 80224121865SAdrian Maldonado 80324121865SAdrian Maldonado PetscFunctionBegin; 80424121865SAdrian Maldonado /* Create index sets to map from "points" to "subpoints" */ 80524121865SAdrian Maldonado ierr = PetscMalloc1(pend - pstart, &subpoints);CHKERRQ(ierr); 80624121865SAdrian Maldonado for (i = pstart; i < pend; i++) { 80724121865SAdrian Maldonado subpoints[i - pstart] = i; 80824121865SAdrian Maldonado } 809459726d8SSatish Balay ierr = ISLocalToGlobalMappingCreate(PETSC_COMM_WORLD,1,pend-pstart,subpoints,PETSC_COPY_VALUES,map);CHKERRQ(ierr); 81024121865SAdrian Maldonado ierr = PetscFree(subpoints);CHKERRQ(ierr); 81124121865SAdrian Maldonado PetscFunctionReturn(0); 81224121865SAdrian Maldonado } 81324121865SAdrian Maldonado 81424121865SAdrian Maldonado /*@ 81524121865SAdrian Maldonado DMNetworkAssembleGraphStructures - Assembles vertex and edge data structures. Must be called after DMNetworkDistribute. 81624121865SAdrian Maldonado 81724121865SAdrian Maldonado Collective 81824121865SAdrian Maldonado 81924121865SAdrian Maldonado Input Parameters: 82024121865SAdrian Maldonado . dm - The DMNetworkObject 82124121865SAdrian Maldonado 82224121865SAdrian Maldonado Note: the routine will create alternative orderings for the vertices and edges. Assume global network points are: 82324121865SAdrian Maldonado 82424121865SAdrian Maldonado points = [0 1 2 3 4 5 6] 82524121865SAdrian Maldonado 82624121865SAdrian Maldonado where edges = [0, 3] and vertices = [4, 6]. The new orderings will be specific to the subset (i.e vertices = [0, 2]). 82724121865SAdrian Maldonado 82824121865SAdrian Maldonado With this new ordering a local PetscSection, global PetscSection and PetscSF will be created specific to the subset. 82924121865SAdrian Maldonado 83024121865SAdrian Maldonado Level: intermediate 83124121865SAdrian Maldonado 83224121865SAdrian Maldonado @*/ 83324121865SAdrian Maldonado PetscErrorCode DMNetworkAssembleGraphStructures(DM dm) 83424121865SAdrian Maldonado { 83524121865SAdrian Maldonado PetscErrorCode ierr; 83624121865SAdrian Maldonado MPI_Comm comm; 8379852e123SBarry Smith PetscMPIInt rank, size; 83824121865SAdrian Maldonado DM_Network *network = (DM_Network*)dm->data; 83924121865SAdrian Maldonado 840eab1376dSHong Zhang PetscFunctionBegin; 84124121865SAdrian Maldonado ierr = PetscObjectGetComm((PetscObject)dm,&comm);CHKERRQ(ierr); 84224121865SAdrian Maldonado ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr); 8439852e123SBarry Smith ierr = MPI_Comm_size(comm, &size);CHKERRQ(ierr); 84424121865SAdrian Maldonado 84524121865SAdrian Maldonado /* Create maps for vertices and edges */ 84624121865SAdrian Maldonado ierr = DMNetworkSetSubMap_private(network->vStart,network->vEnd,&network->vertex.mapping);CHKERRQ(ierr); 84724121865SAdrian Maldonado ierr = DMNetworkSetSubMap_private(network->eStart,network->eEnd,&network->edge.mapping);CHKERRQ(ierr); 84824121865SAdrian Maldonado 84924121865SAdrian Maldonado /* Create local sub-sections */ 85024121865SAdrian Maldonado ierr = DMNetworkGetSubSection_private(network->DofSection,network->vStart,network->vEnd,&network->vertex.DofSection);CHKERRQ(ierr); 85124121865SAdrian Maldonado ierr = DMNetworkGetSubSection_private(network->DofSection,network->eStart,network->eEnd,&network->edge.DofSection);CHKERRQ(ierr); 85224121865SAdrian Maldonado 8539852e123SBarry Smith if (size > 1) { 85424121865SAdrian Maldonado ierr = PetscSFGetSubSF(network->plex->sf, network->vertex.mapping, &network->vertex.sf);CHKERRQ(ierr); 85524121865SAdrian Maldonado ierr = PetscSectionCreateGlobalSection(network->vertex.DofSection, network->vertex.sf, PETSC_FALSE, PETSC_FALSE, &network->vertex.GlobalDofSection);CHKERRQ(ierr); 85624121865SAdrian Maldonado ierr = PetscSFGetSubSF(network->plex->sf, network->edge.mapping, &network->edge.sf);CHKERRQ(ierr); 85724121865SAdrian Maldonado ierr = PetscSectionCreateGlobalSection(network->edge.DofSection, network->edge.sf, PETSC_FALSE, PETSC_FALSE, &network->edge.GlobalDofSection);CHKERRQ(ierr); 85824121865SAdrian Maldonado } else { 85924121865SAdrian Maldonado /* create structures for vertex */ 86024121865SAdrian Maldonado ierr = PetscSectionClone(network->vertex.DofSection,&network->vertex.GlobalDofSection);CHKERRQ(ierr); 86124121865SAdrian Maldonado /* create structures for edge */ 86224121865SAdrian Maldonado ierr = PetscSectionClone(network->edge.DofSection,&network->edge.GlobalDofSection);CHKERRQ(ierr); 86324121865SAdrian Maldonado } 86424121865SAdrian Maldonado 86524121865SAdrian Maldonado 86624121865SAdrian Maldonado /* Add viewers */ 86724121865SAdrian Maldonado ierr = PetscObjectSetName((PetscObject)network->edge.GlobalDofSection,"Global edge dof section");CHKERRQ(ierr); 86824121865SAdrian Maldonado ierr = PetscObjectSetName((PetscObject)network->vertex.GlobalDofSection,"Global vertex dof section");CHKERRQ(ierr); 86924121865SAdrian Maldonado ierr = PetscSectionViewFromOptions(network->edge.GlobalDofSection, NULL, "-edge_global_section_view");CHKERRQ(ierr); 87024121865SAdrian Maldonado ierr = PetscSectionViewFromOptions(network->vertex.GlobalDofSection, NULL, "-vertex_global_section_view");CHKERRQ(ierr); 87124121865SAdrian Maldonado 87224121865SAdrian Maldonado PetscFunctionReturn(0); 87324121865SAdrian Maldonado } 8747b6afd5bSHong Zhang 8755f2c45f1SShri Abhyankar /*@ 8765f2c45f1SShri Abhyankar DMNetworkDistribute - Distributes the network and moves associated component data. 8775f2c45f1SShri Abhyankar 8785f2c45f1SShri Abhyankar Collective 8795f2c45f1SShri Abhyankar 8805f2c45f1SShri Abhyankar Input Parameter: 881d3464fd4SAdrian Maldonado + DM - the DMNetwork object 8825f2c45f1SShri Abhyankar - overlap - The overlap of partitions, 0 is the default 8835f2c45f1SShri Abhyankar 8845f2c45f1SShri Abhyankar Notes: 8858b171c8eSHong Zhang Distributes the network with <overlap>-overlapping partitioning of the edges. 8865f2c45f1SShri Abhyankar 8875f2c45f1SShri Abhyankar Level: intermediate 8885f2c45f1SShri Abhyankar 8895f2c45f1SShri Abhyankar .seealso: DMNetworkCreate 8905f2c45f1SShri Abhyankar @*/ 891d3464fd4SAdrian Maldonado PetscErrorCode DMNetworkDistribute(DM *dm,PetscInt overlap) 8925f2c45f1SShri Abhyankar { 893d3464fd4SAdrian Maldonado MPI_Comm comm; 8945f2c45f1SShri Abhyankar PetscErrorCode ierr; 895d3464fd4SAdrian Maldonado PetscMPIInt size; 896d3464fd4SAdrian Maldonado DM_Network *oldDMnetwork = (DM_Network*)((*dm)->data); 897d3464fd4SAdrian Maldonado DM_Network *newDMnetwork; 8985f2c45f1SShri Abhyankar PetscSF pointsf; 8995f2c45f1SShri Abhyankar DM newDM; 90051ac5effSHong Zhang PetscPartitioner part; 9015f2c45f1SShri Abhyankar 9025f2c45f1SShri Abhyankar PetscFunctionBegin; 903d3464fd4SAdrian Maldonado 904d3464fd4SAdrian Maldonado ierr = PetscObjectGetComm((PetscObject)*dm,&comm);CHKERRQ(ierr); 905d3464fd4SAdrian Maldonado ierr = MPI_Comm_size(comm, &size);CHKERRQ(ierr); 906d3464fd4SAdrian Maldonado if (size == 1) PetscFunctionReturn(0); 907d3464fd4SAdrian Maldonado 908d3464fd4SAdrian Maldonado ierr = DMNetworkCreate(PetscObjectComm((PetscObject)*dm),&newDM);CHKERRQ(ierr); 9095f2c45f1SShri Abhyankar newDMnetwork = (DM_Network*)newDM->data; 9105f2c45f1SShri Abhyankar newDMnetwork->dataheadersize = sizeof(struct _p_DMNetworkComponentHeader)/sizeof(DMNetworkComponentGenericDataType); 91151ac5effSHong Zhang 91251ac5effSHong Zhang /* Enable runtime options for petscpartitioner */ 91351ac5effSHong Zhang ierr = DMPlexGetPartitioner(oldDMnetwork->plex,&part);CHKERRQ(ierr); 91451ac5effSHong Zhang ierr = PetscPartitionerSetFromOptions(part);CHKERRQ(ierr); 91551ac5effSHong Zhang 9165f2c45f1SShri Abhyankar /* Distribute plex dm and dof section */ 91780cf41d5SMatthew G. Knepley ierr = DMPlexDistribute(oldDMnetwork->plex,overlap,&pointsf,&newDMnetwork->plex);CHKERRQ(ierr); 91851ac5effSHong Zhang 9195f2c45f1SShri Abhyankar /* Distribute dof section */ 920d3464fd4SAdrian Maldonado ierr = PetscSectionCreate(PetscObjectComm((PetscObject)*dm),&newDMnetwork->DofSection);CHKERRQ(ierr); 9215f2c45f1SShri Abhyankar ierr = PetscSFDistributeSection(pointsf,oldDMnetwork->DofSection,NULL,newDMnetwork->DofSection);CHKERRQ(ierr); 922d3464fd4SAdrian Maldonado ierr = PetscSectionCreate(PetscObjectComm((PetscObject)*dm),&newDMnetwork->DataSection);CHKERRQ(ierr); 92351ac5effSHong Zhang 9245f2c45f1SShri Abhyankar /* Distribute data and associated section */ 92531da1fc8SHong Zhang ierr = DMPlexDistributeData(newDMnetwork->plex,pointsf,oldDMnetwork->DataSection,MPIU_INT,(void*)oldDMnetwork->componentdataarray,newDMnetwork->DataSection,(void**)&newDMnetwork->componentdataarray);CHKERRQ(ierr); 92624121865SAdrian Maldonado 9275f2c45f1SShri Abhyankar ierr = PetscSectionGetChart(newDMnetwork->DataSection,&newDMnetwork->pStart,&newDMnetwork->pEnd);CHKERRQ(ierr); 9285f2c45f1SShri Abhyankar ierr = DMPlexGetHeightStratum(newDMnetwork->plex,0, &newDMnetwork->eStart,&newDMnetwork->eEnd);CHKERRQ(ierr); 9295f2c45f1SShri Abhyankar ierr = DMPlexGetHeightStratum(newDMnetwork->plex,1,&newDMnetwork->vStart,&newDMnetwork->vEnd);CHKERRQ(ierr); 9305f2c45f1SShri Abhyankar newDMnetwork->nEdges = newDMnetwork->eEnd - newDMnetwork->eStart; 9316fefedf4SHong Zhang newDMnetwork->nVertices = newDMnetwork->vEnd - newDMnetwork->vStart; 9326fefedf4SHong Zhang newDMnetwork->NVertices = oldDMnetwork->NVertices; 9335f2c45f1SShri Abhyankar newDMnetwork->NEdges = oldDMnetwork->NEdges; 93424121865SAdrian Maldonado 9355f2c45f1SShri Abhyankar /* Set Dof section as the default section for dm */ 9365f2c45f1SShri Abhyankar ierr = DMSetDefaultSection(newDMnetwork->plex,newDMnetwork->DofSection);CHKERRQ(ierr); 9375f2c45f1SShri Abhyankar ierr = DMGetDefaultGlobalSection(newDMnetwork->plex,&newDMnetwork->GlobalDofSection);CHKERRQ(ierr); 9385f2c45f1SShri Abhyankar 93924121865SAdrian Maldonado /* Destroy point SF */ 94024121865SAdrian Maldonado ierr = PetscSFDestroy(&pointsf);CHKERRQ(ierr); 94124121865SAdrian Maldonado 942d3464fd4SAdrian Maldonado ierr = DMDestroy(dm);CHKERRQ(ierr); 943d3464fd4SAdrian Maldonado *dm = newDM; 9445f2c45f1SShri Abhyankar PetscFunctionReturn(0); 9455f2c45f1SShri Abhyankar } 9465f2c45f1SShri Abhyankar 94724121865SAdrian Maldonado /*@C 94824121865SAdrian Maldonado PetscSFGetSubSF - Returns an SF for a specific subset of points. Leaves are re-numbered to reflect the new ordering. 94924121865SAdrian Maldonado 95024121865SAdrian Maldonado Input Parameters: 95124121865SAdrian Maldonado + masterSF - the original SF structure 95224121865SAdrian Maldonado - map - a ISLocalToGlobal mapping that contains the subset of points 95324121865SAdrian Maldonado 95424121865SAdrian Maldonado Output Parameters: 95524121865SAdrian Maldonado . subSF - a subset of the masterSF for the desired subset. 95624121865SAdrian Maldonado */ 95724121865SAdrian Maldonado 95824121865SAdrian Maldonado PetscErrorCode PetscSFGetSubSF(PetscSF mastersf, ISLocalToGlobalMapping map, PetscSF *subSF) { 95924121865SAdrian Maldonado 96024121865SAdrian Maldonado PetscErrorCode ierr; 96124121865SAdrian Maldonado PetscInt nroots, nleaves, *ilocal_sub; 96224121865SAdrian Maldonado PetscInt i, *ilocal_map, nroots_sub, nleaves_sub = 0; 96324121865SAdrian Maldonado PetscInt *local_points, *remote_points; 96424121865SAdrian Maldonado PetscSFNode *iremote_sub; 96524121865SAdrian Maldonado const PetscInt *ilocal; 96624121865SAdrian Maldonado const PetscSFNode *iremote; 96724121865SAdrian Maldonado 96824121865SAdrian Maldonado PetscFunctionBegin; 96924121865SAdrian Maldonado ierr = PetscSFGetGraph(mastersf,&nroots,&nleaves,&ilocal,&iremote);CHKERRQ(ierr); 97024121865SAdrian Maldonado 97124121865SAdrian Maldonado /* Look for leaves that pertain to the subset of points. Get the local ordering */ 97224121865SAdrian Maldonado ierr = PetscMalloc1(nleaves,&ilocal_map);CHKERRQ(ierr); 97324121865SAdrian Maldonado ierr = ISGlobalToLocalMappingApply(map,IS_GTOLM_MASK,nleaves,ilocal,NULL,ilocal_map);CHKERRQ(ierr); 97424121865SAdrian Maldonado for (i = 0; i < nleaves; i++) { 97524121865SAdrian Maldonado if (ilocal_map[i] != -1) nleaves_sub += 1; 97624121865SAdrian Maldonado } 97724121865SAdrian Maldonado /* Re-number ilocal with subset numbering. Need information from roots */ 97824121865SAdrian Maldonado ierr = PetscMalloc2(nroots,&local_points,nroots,&remote_points);CHKERRQ(ierr); 97924121865SAdrian Maldonado for (i = 0; i < nroots; i++) local_points[i] = i; 98024121865SAdrian Maldonado ierr = ISGlobalToLocalMappingApply(map,IS_GTOLM_MASK,nroots,local_points,NULL,local_points);CHKERRQ(ierr); 98124121865SAdrian Maldonado ierr = PetscSFBcastBegin(mastersf, MPIU_INT, local_points, remote_points);CHKERRQ(ierr); 98224121865SAdrian Maldonado ierr = PetscSFBcastEnd(mastersf, MPIU_INT, local_points, remote_points);CHKERRQ(ierr); 98324121865SAdrian Maldonado /* Fill up graph using local (that is, local to the subset) numbering. */ 9844b70a8deSAdrian Maldonado ierr = PetscMalloc1(nleaves_sub,&ilocal_sub);CHKERRQ(ierr); 9854b70a8deSAdrian Maldonado ierr = PetscMalloc1(nleaves_sub,&iremote_sub);CHKERRQ(ierr); 98624121865SAdrian Maldonado nleaves_sub = 0; 98724121865SAdrian Maldonado for (i = 0; i < nleaves; i++) { 98824121865SAdrian Maldonado if (ilocal_map[i] != -1) { 98924121865SAdrian Maldonado ilocal_sub[nleaves_sub] = ilocal_map[i]; 9904b70a8deSAdrian Maldonado iremote_sub[nleaves_sub].rank = iremote[i].rank; 99124121865SAdrian Maldonado iremote_sub[nleaves_sub].index = remote_points[ilocal[i]]; 99224121865SAdrian Maldonado nleaves_sub += 1; 99324121865SAdrian Maldonado } 99424121865SAdrian Maldonado } 99524121865SAdrian Maldonado ierr = PetscFree2(local_points,remote_points);CHKERRQ(ierr); 99624121865SAdrian Maldonado ierr = ISLocalToGlobalMappingGetSize(map,&nroots_sub);CHKERRQ(ierr); 99724121865SAdrian Maldonado 99824121865SAdrian Maldonado /* Create new subSF */ 99924121865SAdrian Maldonado ierr = PetscSFCreate(PETSC_COMM_WORLD,subSF);CHKERRQ(ierr); 100024121865SAdrian Maldonado ierr = PetscSFSetFromOptions(*subSF);CHKERRQ(ierr); 10014b70a8deSAdrian Maldonado ierr = PetscSFSetGraph(*subSF,nroots_sub,nleaves_sub,ilocal_sub,PETSC_OWN_POINTER,iremote_sub,PETSC_COPY_VALUES);CHKERRQ(ierr); 100224121865SAdrian Maldonado ierr = PetscFree(ilocal_map);CHKERRQ(ierr); 10034b70a8deSAdrian Maldonado ierr = PetscFree(iremote_sub);CHKERRQ(ierr); 100424121865SAdrian Maldonado PetscFunctionReturn(0); 100524121865SAdrian Maldonado } 100624121865SAdrian Maldonado 10075f2c45f1SShri Abhyankar /*@C 10085f2c45f1SShri Abhyankar DMNetworkGetSupportingEdges - Return the supporting edges for this vertex point 10095f2c45f1SShri Abhyankar 10105f2c45f1SShri Abhyankar Not Collective 10115f2c45f1SShri Abhyankar 10125f2c45f1SShri Abhyankar Input Parameters: 10135f2c45f1SShri Abhyankar + dm - The DMNetwork object 10145f2c45f1SShri Abhyankar - p - the vertex point 10155f2c45f1SShri Abhyankar 10165f2c45f1SShri Abhyankar Output Paramters: 10175f2c45f1SShri Abhyankar + nedges - number of edges connected to this vertex point 10185f2c45f1SShri Abhyankar - edges - List of edge points 10195f2c45f1SShri Abhyankar 10205f2c45f1SShri Abhyankar Level: intermediate 10215f2c45f1SShri Abhyankar 10225f2c45f1SShri Abhyankar Fortran Notes: 10235f2c45f1SShri Abhyankar Since it returns an array, this routine is only available in Fortran 90, and you must 10245f2c45f1SShri Abhyankar include petsc.h90 in your code. 10255f2c45f1SShri Abhyankar 1026d842c372SHong Zhang .seealso: DMNetworkCreate, DMNetworkGetConnectedVertices 10275f2c45f1SShri Abhyankar @*/ 10285f2c45f1SShri Abhyankar PetscErrorCode DMNetworkGetSupportingEdges(DM dm,PetscInt vertex,PetscInt *nedges,const PetscInt *edges[]) 10295f2c45f1SShri Abhyankar { 10305f2c45f1SShri Abhyankar PetscErrorCode ierr; 10315f2c45f1SShri Abhyankar DM_Network *network = (DM_Network*)dm->data; 10325f2c45f1SShri Abhyankar 10335f2c45f1SShri Abhyankar PetscFunctionBegin; 10345f2c45f1SShri Abhyankar ierr = DMPlexGetSupportSize(network->plex,vertex,nedges);CHKERRQ(ierr); 10355f2c45f1SShri Abhyankar ierr = DMPlexGetSupport(network->plex,vertex,edges);CHKERRQ(ierr); 10365f2c45f1SShri Abhyankar PetscFunctionReturn(0); 10375f2c45f1SShri Abhyankar } 10385f2c45f1SShri Abhyankar 10395f2c45f1SShri Abhyankar /*@C 1040d842c372SHong Zhang DMNetworkGetConnectedVertices - Return the connected vertices for this edge point 10415f2c45f1SShri Abhyankar 10425f2c45f1SShri Abhyankar Not Collective 10435f2c45f1SShri Abhyankar 10445f2c45f1SShri Abhyankar Input Parameters: 10455f2c45f1SShri Abhyankar + dm - The DMNetwork object 10465f2c45f1SShri Abhyankar - p - the edge point 10475f2c45f1SShri Abhyankar 10485f2c45f1SShri Abhyankar Output Paramters: 10495f2c45f1SShri Abhyankar . vertices - vertices connected to this edge 10505f2c45f1SShri Abhyankar 10515f2c45f1SShri Abhyankar Level: intermediate 10525f2c45f1SShri Abhyankar 10535f2c45f1SShri Abhyankar Fortran Notes: 10545f2c45f1SShri Abhyankar Since it returns an array, this routine is only available in Fortran 90, and you must 10555f2c45f1SShri Abhyankar include petsc.h90 in your code. 10565f2c45f1SShri Abhyankar 10575f2c45f1SShri Abhyankar .seealso: DMNetworkCreate, DMNetworkGetSupportingEdges 10585f2c45f1SShri Abhyankar @*/ 1059d842c372SHong Zhang PetscErrorCode DMNetworkGetConnectedVertices(DM dm,PetscInt edge,const PetscInt *vertices[]) 10605f2c45f1SShri Abhyankar { 10615f2c45f1SShri Abhyankar PetscErrorCode ierr; 10625f2c45f1SShri Abhyankar DM_Network *network = (DM_Network*)dm->data; 10635f2c45f1SShri Abhyankar 10645f2c45f1SShri Abhyankar PetscFunctionBegin; 10655f2c45f1SShri Abhyankar ierr = DMPlexGetCone(network->plex,edge,vertices);CHKERRQ(ierr); 10665f2c45f1SShri Abhyankar PetscFunctionReturn(0); 10675f2c45f1SShri Abhyankar } 10685f2c45f1SShri Abhyankar 10695f2c45f1SShri Abhyankar /*@ 10705f2c45f1SShri Abhyankar DMNetworkIsGhostVertex - Returns TRUE if the vertex is a ghost vertex 10715f2c45f1SShri Abhyankar 10725f2c45f1SShri Abhyankar Not Collective 10735f2c45f1SShri Abhyankar 10745f2c45f1SShri Abhyankar Input Parameters: 10755f2c45f1SShri Abhyankar + dm - The DMNetwork object 10765f2c45f1SShri Abhyankar . p - the vertex point 10775f2c45f1SShri Abhyankar 10785f2c45f1SShri Abhyankar Output Parameter: 10795f2c45f1SShri Abhyankar . isghost - TRUE if the vertex is a ghost point 10805f2c45f1SShri Abhyankar 10815f2c45f1SShri Abhyankar Level: intermediate 10825f2c45f1SShri Abhyankar 1083d842c372SHong Zhang .seealso: DMNetworkCreate, DMNetworkGetConnectedVertices, DMNetworkGetVertexRange 10845f2c45f1SShri Abhyankar @*/ 10855f2c45f1SShri Abhyankar PetscErrorCode DMNetworkIsGhostVertex(DM dm,PetscInt p,PetscBool *isghost) 10865f2c45f1SShri Abhyankar { 10875f2c45f1SShri Abhyankar PetscErrorCode ierr; 10885f2c45f1SShri Abhyankar DM_Network *network = (DM_Network*)dm->data; 10895f2c45f1SShri Abhyankar PetscInt offsetg; 10905f2c45f1SShri Abhyankar PetscSection sectiong; 10915f2c45f1SShri Abhyankar 10925f2c45f1SShri Abhyankar PetscFunctionBegin; 10935f2c45f1SShri Abhyankar *isghost = PETSC_FALSE; 10945f2c45f1SShri Abhyankar ierr = DMGetDefaultGlobalSection(network->plex,§iong);CHKERRQ(ierr); 10955f2c45f1SShri Abhyankar ierr = PetscSectionGetOffset(sectiong,p,&offsetg);CHKERRQ(ierr); 10965f2c45f1SShri Abhyankar if (offsetg < 0) *isghost = PETSC_TRUE; 10975f2c45f1SShri Abhyankar PetscFunctionReturn(0); 10985f2c45f1SShri Abhyankar } 10995f2c45f1SShri Abhyankar 11005f2c45f1SShri Abhyankar PetscErrorCode DMSetUp_Network(DM dm) 11015f2c45f1SShri Abhyankar { 11025f2c45f1SShri Abhyankar PetscErrorCode ierr; 11035f2c45f1SShri Abhyankar DM_Network *network=(DM_Network*)dm->data; 11045f2c45f1SShri Abhyankar 11055f2c45f1SShri Abhyankar PetscFunctionBegin; 11065f2c45f1SShri Abhyankar ierr = DMNetworkComponentSetUp(dm);CHKERRQ(ierr); 11075f2c45f1SShri Abhyankar ierr = DMNetworkVariablesSetUp(dm);CHKERRQ(ierr); 11085f2c45f1SShri Abhyankar 11095f2c45f1SShri Abhyankar ierr = DMSetDefaultSection(network->plex,network->DofSection);CHKERRQ(ierr); 11105f2c45f1SShri Abhyankar ierr = DMGetDefaultGlobalSection(network->plex,&network->GlobalDofSection);CHKERRQ(ierr); 11115f2c45f1SShri Abhyankar PetscFunctionReturn(0); 11125f2c45f1SShri Abhyankar } 11135f2c45f1SShri Abhyankar 11141ad426b7SHong Zhang /*@ 111517df6e9eSHong Zhang DMNetworkHasJacobian - Sets global flag for using user's sub Jacobian matrices 11161ad426b7SHong Zhang -- replaced by DMNetworkSetOption(network,userjacobian,PETSC_TURE)? 11171ad426b7SHong Zhang 11181ad426b7SHong Zhang Collective 11191ad426b7SHong Zhang 11201ad426b7SHong Zhang Input Parameters: 112183b2e829SHong Zhang + dm - The DMNetwork object 112283b2e829SHong Zhang . eflg - turn the option on (PETSC_TRUE) or off (PETSC_FALSE) if user provides Jacobian for edges 112383b2e829SHong Zhang - vflg - turn the option on (PETSC_TRUE) or off (PETSC_FALSE) if user provides Jacobian for vertices 11241ad426b7SHong Zhang 11251ad426b7SHong Zhang Level: intermediate 11261ad426b7SHong Zhang 11271ad426b7SHong Zhang @*/ 112883b2e829SHong Zhang PetscErrorCode DMNetworkHasJacobian(DM dm,PetscBool eflg,PetscBool vflg) 11291ad426b7SHong Zhang { 11301ad426b7SHong Zhang DM_Network *network=(DM_Network*)dm->data; 11311ad426b7SHong Zhang 11321ad426b7SHong Zhang PetscFunctionBegin; 113383b2e829SHong Zhang network->userEdgeJacobian = eflg; 113483b2e829SHong Zhang network->userVertexJacobian = vflg; 11351ad426b7SHong Zhang PetscFunctionReturn(0); 11361ad426b7SHong Zhang } 11371ad426b7SHong Zhang 11381ad426b7SHong Zhang /*@ 113983b2e829SHong Zhang DMNetworkEdgeSetMatrix - Sets user-provided Jacobian matrices for this edge to the network 114083b2e829SHong Zhang 114183b2e829SHong Zhang Not Collective 114283b2e829SHong Zhang 114383b2e829SHong Zhang Input Parameters: 114483b2e829SHong Zhang + dm - The DMNetwork object 114583b2e829SHong Zhang . p - the edge point 11463e97b6e8SHong Zhang - J - array (size = 3) of Jacobian submatrices for this edge point: 11473e97b6e8SHong Zhang J[0]: this edge 1148d842c372SHong Zhang J[1] and J[2]: connected vertices, obtained by calling DMNetworkGetConnectedVertices() 114983b2e829SHong Zhang 115083b2e829SHong Zhang Level: intermediate 115183b2e829SHong Zhang 115283b2e829SHong Zhang .seealso: DMNetworkVertexSetMatrix 115383b2e829SHong Zhang @*/ 115483b2e829SHong Zhang PetscErrorCode DMNetworkEdgeSetMatrix(DM dm,PetscInt p,Mat J[]) 115583b2e829SHong Zhang { 115683b2e829SHong Zhang PetscErrorCode ierr; 115783b2e829SHong Zhang DM_Network *network=(DM_Network*)dm->data; 115883b2e829SHong Zhang 115983b2e829SHong Zhang PetscFunctionBegin; 1160883e35e8SHong Zhang if (!network->userEdgeJacobian) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ORDER,"Must call DMNetworkHasJacobian() collectively before calling DMNetworkEdgeSetMatrix"); 116183b2e829SHong Zhang if (!network->Je) { 1162883e35e8SHong Zhang ierr = PetscCalloc1(3*network->nEdges,&network->Je);CHKERRQ(ierr); 116383b2e829SHong Zhang } 1164883e35e8SHong Zhang network->Je[3*p] = J[0]; 1165883e35e8SHong Zhang network->Je[3*p+1] = J[1]; 1166883e35e8SHong Zhang network->Je[3*p+2] = J[2]; 116783b2e829SHong Zhang PetscFunctionReturn(0); 116883b2e829SHong Zhang } 116983b2e829SHong Zhang 117083b2e829SHong Zhang /*@ 117176ddfea5SHong Zhang DMNetworkVertexSetMatrix - Sets user-provided Jacobian matrix for this vertex to the network 11721ad426b7SHong Zhang 11731ad426b7SHong Zhang Not Collective 11741ad426b7SHong Zhang 11751ad426b7SHong Zhang Input Parameters: 11761ad426b7SHong Zhang + dm - The DMNetwork object 11771ad426b7SHong Zhang . p - the vertex point 11783e97b6e8SHong Zhang - J - array of Jacobian (size = 2*(num of supporting edges) + 1) submatrices for this vertex point: 11793e97b6e8SHong Zhang J[0]: this vertex 11803e97b6e8SHong Zhang J[1+2*i]: i-th supporting edge 11813e97b6e8SHong Zhang J[1+2*i+1]: i-th connected vertex 11821ad426b7SHong Zhang 11831ad426b7SHong Zhang Level: intermediate 11841ad426b7SHong Zhang 118583b2e829SHong Zhang .seealso: DMNetworkEdgeSetMatrix 11861ad426b7SHong Zhang @*/ 1187883e35e8SHong Zhang PetscErrorCode DMNetworkVertexSetMatrix(DM dm,PetscInt p,Mat J[]) 11885f2c45f1SShri Abhyankar { 11895f2c45f1SShri Abhyankar PetscErrorCode ierr; 11905f2c45f1SShri Abhyankar DM_Network *network=(DM_Network*)dm->data; 1191f4431b8cSHong Zhang PetscInt i,*vptr,nedges,vStart,vEnd; 1192883e35e8SHong Zhang const PetscInt *edges; 11935f2c45f1SShri Abhyankar 11945f2c45f1SShri Abhyankar PetscFunctionBegin; 1195883e35e8SHong Zhang if (!network->userVertexJacobian) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ORDER,"Must call DMNetworkHasJacobian() collectively before calling DMNetworkVertexSetMatrix"); 1196883e35e8SHong Zhang 1197883e35e8SHong Zhang ierr = DMNetworkGetVertexRange(dm,&vStart,&vEnd);CHKERRQ(ierr); 1198f4431b8cSHong Zhang 119983b2e829SHong Zhang if (!network->Jv) { 12006fefedf4SHong Zhang PetscInt nVertices = network->nVertices,nedges_total; 1201883e35e8SHong Zhang /* count nvertex_total */ 12023e97b6e8SHong Zhang nedges_total = 0; 1203883e35e8SHong Zhang ierr = DMNetworkGetVertexRange(dm,&vStart,&vEnd);CHKERRQ(ierr); 12046fefedf4SHong Zhang ierr = PetscMalloc1(nVertices+1,&vptr);CHKERRQ(ierr); 1205f4431b8cSHong Zhang 1206883e35e8SHong Zhang vptr[0] = 0; 12076fefedf4SHong Zhang for (i=0; i<nVertices; i++) { 1208f4431b8cSHong Zhang ierr = DMNetworkGetSupportingEdges(dm,i+vStart,&nedges,&edges);CHKERRQ(ierr); 1209883e35e8SHong Zhang nedges_total += nedges; 1210f4431b8cSHong Zhang vptr[i+1] = vptr[i] + 2*nedges + 1; 12111ad426b7SHong Zhang } 12123e97b6e8SHong Zhang 12136fefedf4SHong Zhang ierr = PetscCalloc1(2*nedges_total+nVertices,&network->Jv);CHKERRQ(ierr); 1214883e35e8SHong Zhang network->Jvptr = vptr; 1215883e35e8SHong Zhang } 1216883e35e8SHong Zhang 1217883e35e8SHong Zhang vptr = network->Jvptr; 12183e97b6e8SHong Zhang network->Jv[vptr[p-vStart]] = J[0]; /* Set Jacobian for this vertex */ 12193e97b6e8SHong Zhang 12203e97b6e8SHong Zhang /* Set Jacobian for each supporting edge and connected vertex */ 1221883e35e8SHong Zhang ierr = DMNetworkGetSupportingEdges(dm,p,&nedges,&edges);CHKERRQ(ierr); 1222883e35e8SHong Zhang for (i=1; i<=2*nedges; i++) network->Jv[vptr[p-vStart]+i] = J[i]; 1223883e35e8SHong Zhang PetscFunctionReturn(0); 1224883e35e8SHong Zhang } 1225883e35e8SHong Zhang 1226e0f69777SHong Zhang PETSC_STATIC_INLINE PetscErrorCode MatSetPreallocationDenseblock_private(PetscInt nrows,PetscInt *rows,PetscInt ncols,PetscBool ghost,Vec vdnz,Vec vonz) 12275cf7da58SHong Zhang { 12285cf7da58SHong Zhang PetscErrorCode ierr; 12295cf7da58SHong Zhang PetscInt j; 12305cf7da58SHong Zhang PetscScalar val=(PetscScalar)ncols; 12315cf7da58SHong Zhang 12325cf7da58SHong Zhang PetscFunctionBegin; 12335cf7da58SHong Zhang if (!ghost) { 12345cf7da58SHong Zhang for (j=0; j<nrows; j++) { 12355cf7da58SHong Zhang ierr = VecSetValues(vdnz,1,&rows[j],&val,ADD_VALUES);CHKERRQ(ierr); 12365cf7da58SHong Zhang } 12375cf7da58SHong Zhang } else { 12385cf7da58SHong Zhang for (j=0; j<nrows; j++) { 12395cf7da58SHong Zhang ierr = VecSetValues(vonz,1,&rows[j],&val,ADD_VALUES);CHKERRQ(ierr); 12405cf7da58SHong Zhang } 12415cf7da58SHong Zhang } 12425cf7da58SHong Zhang PetscFunctionReturn(0); 12435cf7da58SHong Zhang } 12445cf7da58SHong Zhang 1245e0f69777SHong Zhang PETSC_STATIC_INLINE PetscErrorCode MatSetPreallocationUserblock_private(Mat Ju,PetscInt nrows,PetscInt *rows,PetscInt ncols,PetscBool ghost,Vec vdnz,Vec vonz) 12465cf7da58SHong Zhang { 12475cf7da58SHong Zhang PetscErrorCode ierr; 12485cf7da58SHong Zhang PetscInt j,ncols_u; 12495cf7da58SHong Zhang PetscScalar val; 12505cf7da58SHong Zhang 12515cf7da58SHong Zhang PetscFunctionBegin; 12525cf7da58SHong Zhang if (!ghost) { 12535cf7da58SHong Zhang for (j=0; j<nrows; j++) { 12545cf7da58SHong Zhang ierr = MatGetRow(Ju,j,&ncols_u,NULL,NULL);CHKERRQ(ierr); 12555cf7da58SHong Zhang val = (PetscScalar)ncols_u; 12565cf7da58SHong Zhang ierr = VecSetValues(vdnz,1,&rows[j],&val,ADD_VALUES);CHKERRQ(ierr); 12575cf7da58SHong Zhang ierr = MatRestoreRow(Ju,j,&ncols_u,NULL,NULL);CHKERRQ(ierr); 12585cf7da58SHong Zhang } 12595cf7da58SHong Zhang } else { 12605cf7da58SHong Zhang for (j=0; j<nrows; j++) { 12615cf7da58SHong Zhang ierr = MatGetRow(Ju,j,&ncols_u,NULL,NULL);CHKERRQ(ierr); 12625cf7da58SHong Zhang val = (PetscScalar)ncols_u; 12635cf7da58SHong Zhang ierr = VecSetValues(vonz,1,&rows[j],&val,ADD_VALUES);CHKERRQ(ierr); 12645cf7da58SHong Zhang ierr = MatRestoreRow(Ju,j,&ncols_u,NULL,NULL);CHKERRQ(ierr); 12655cf7da58SHong Zhang } 12665cf7da58SHong Zhang } 12675cf7da58SHong Zhang PetscFunctionReturn(0); 12685cf7da58SHong Zhang } 12695cf7da58SHong Zhang 1270e0f69777SHong Zhang PETSC_STATIC_INLINE PetscErrorCode MatSetPreallocationblock_private(Mat Ju,PetscInt nrows,PetscInt *rows,PetscInt ncols,PetscBool ghost,Vec vdnz,Vec vonz) 12715cf7da58SHong Zhang { 12725cf7da58SHong Zhang PetscErrorCode ierr; 12735cf7da58SHong Zhang 12745cf7da58SHong Zhang PetscFunctionBegin; 12755cf7da58SHong Zhang if (Ju) { 12765cf7da58SHong Zhang ierr = MatSetPreallocationUserblock_private(Ju,nrows,rows,ncols,ghost,vdnz,vonz);CHKERRQ(ierr); 12775cf7da58SHong Zhang } else { 12785cf7da58SHong Zhang ierr = MatSetPreallocationDenseblock_private(nrows,rows,ncols,ghost,vdnz,vonz);CHKERRQ(ierr); 12795cf7da58SHong Zhang } 12805cf7da58SHong Zhang PetscFunctionReturn(0); 12815cf7da58SHong Zhang } 12825cf7da58SHong Zhang 1283e0f69777SHong Zhang PETSC_STATIC_INLINE PetscErrorCode MatSetDenseblock_private(PetscInt nrows,PetscInt *rows,PetscInt ncols,PetscInt cstart,Mat *J) 1284883e35e8SHong Zhang { 1285883e35e8SHong Zhang PetscErrorCode ierr; 1286883e35e8SHong Zhang PetscInt j,*cols; 1287883e35e8SHong Zhang PetscScalar *zeros; 1288883e35e8SHong Zhang 1289883e35e8SHong Zhang PetscFunctionBegin; 1290883e35e8SHong Zhang ierr = PetscCalloc2(ncols,&cols,nrows*ncols,&zeros);CHKERRQ(ierr); 1291883e35e8SHong Zhang for (j=0; j<ncols; j++) cols[j] = j+ cstart; 1292883e35e8SHong Zhang ierr = MatSetValues(*J,nrows,rows,ncols,cols,zeros,INSERT_VALUES);CHKERRQ(ierr); 1293883e35e8SHong Zhang ierr = PetscFree2(cols,zeros);CHKERRQ(ierr); 12941ad426b7SHong Zhang PetscFunctionReturn(0); 12951ad426b7SHong Zhang } 1296a4e85ca8SHong Zhang 1297e0f69777SHong Zhang PETSC_STATIC_INLINE PetscErrorCode MatSetUserblock_private(Mat Ju,PetscInt nrows,PetscInt *rows,PetscInt ncols,PetscInt cstart,Mat *J) 12983e97b6e8SHong Zhang { 12993e97b6e8SHong Zhang PetscErrorCode ierr; 13003e97b6e8SHong Zhang PetscInt j,M,N,row,col,ncols_u; 13013e97b6e8SHong Zhang const PetscInt *cols; 13023e97b6e8SHong Zhang PetscScalar zero=0.0; 13033e97b6e8SHong Zhang 13043e97b6e8SHong Zhang PetscFunctionBegin; 13053e97b6e8SHong Zhang ierr = MatGetSize(Ju,&M,&N);CHKERRQ(ierr); 13063e97b6e8SHong Zhang if (nrows != M || ncols != N) SETERRQ4(PetscObjectComm((PetscObject)Ju),PETSC_ERR_USER,"%D by %D must equal %D by %D",nrows,ncols,M,N); 13073e97b6e8SHong Zhang 13083e97b6e8SHong Zhang for (row=0; row<nrows; row++) { 13093e97b6e8SHong Zhang ierr = MatGetRow(Ju,row,&ncols_u,&cols,NULL);CHKERRQ(ierr); 13103e97b6e8SHong Zhang for (j=0; j<ncols_u; j++) { 13113e97b6e8SHong Zhang col = cols[j] + cstart; 13123e97b6e8SHong Zhang ierr = MatSetValues(*J,1,&rows[row],1,&col,&zero,INSERT_VALUES);CHKERRQ(ierr); 13133e97b6e8SHong Zhang } 13143e97b6e8SHong Zhang ierr = MatRestoreRow(Ju,row,&ncols_u,&cols,NULL);CHKERRQ(ierr); 13153e97b6e8SHong Zhang } 13163e97b6e8SHong Zhang PetscFunctionReturn(0); 13173e97b6e8SHong Zhang } 13181ad426b7SHong Zhang 1319e0f69777SHong Zhang PETSC_STATIC_INLINE PetscErrorCode MatSetblock_private(Mat Ju,PetscInt nrows,PetscInt *rows,PetscInt ncols,PetscInt cstart,Mat *J) 1320a4e85ca8SHong Zhang { 1321a4e85ca8SHong Zhang PetscErrorCode ierr; 1322f4431b8cSHong Zhang 1323a4e85ca8SHong Zhang PetscFunctionBegin; 1324a4e85ca8SHong Zhang if (Ju) { 1325a4e85ca8SHong Zhang ierr = MatSetUserblock_private(Ju,nrows,rows,ncols,cstart,J);CHKERRQ(ierr); 1326a4e85ca8SHong Zhang } else { 1327a4e85ca8SHong Zhang ierr = MatSetDenseblock_private(nrows,rows,ncols,cstart,J);CHKERRQ(ierr); 1328a4e85ca8SHong Zhang } 1329a4e85ca8SHong Zhang PetscFunctionReturn(0); 1330a4e85ca8SHong Zhang } 1331a4e85ca8SHong Zhang 133224121865SAdrian Maldonado /* Creates a GlobalToLocal mapping with a Local and Global section. This is akin to the routine DMGetLocalToGlobalMapping but without the need of providing a dm. 133324121865SAdrian Maldonado */ 133424121865SAdrian Maldonado PetscErrorCode CreateSubGlobalToLocalMapping_private(PetscSection globalsec, PetscSection localsec, ISLocalToGlobalMapping *ltog) 133524121865SAdrian Maldonado { 133624121865SAdrian Maldonado PetscErrorCode ierr; 133724121865SAdrian Maldonado PetscInt i, size, dof; 133824121865SAdrian Maldonado PetscInt *glob2loc; 133924121865SAdrian Maldonado 134024121865SAdrian Maldonado PetscFunctionBegin; 134124121865SAdrian Maldonado ierr = PetscSectionGetStorageSize(localsec,&size);CHKERRQ(ierr); 134224121865SAdrian Maldonado ierr = PetscMalloc1(size,&glob2loc);CHKERRQ(ierr); 134324121865SAdrian Maldonado 134424121865SAdrian Maldonado for (i = 0; i < size; i++) { 134524121865SAdrian Maldonado ierr = PetscSectionGetOffset(globalsec,i,&dof);CHKERRQ(ierr); 134624121865SAdrian Maldonado dof = (dof >= 0) ? dof : -(dof + 1); 134724121865SAdrian Maldonado glob2loc[i] = dof; 134824121865SAdrian Maldonado } 134924121865SAdrian Maldonado 135024121865SAdrian Maldonado ierr = ISLocalToGlobalMappingCreate(PETSC_COMM_WORLD,1,size,glob2loc,PETSC_OWN_POINTER,ltog);CHKERRQ(ierr); 135124121865SAdrian Maldonado #if 0 135224121865SAdrian Maldonado ierr = PetscIntView(size,glob2loc,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 135324121865SAdrian Maldonado #endif 135424121865SAdrian Maldonado PetscFunctionReturn(0); 135524121865SAdrian Maldonado } 135624121865SAdrian Maldonado 135701ad2aeeSHong Zhang #include <petsc/private/matimpl.h> 13581ad426b7SHong Zhang PetscErrorCode DMCreateMatrix_Network(DM dm,Mat *J) 13591ad426b7SHong Zhang { 13601ad426b7SHong Zhang PetscErrorCode ierr; 136124121865SAdrian Maldonado PetscMPIInt rank, size; 13621ad426b7SHong Zhang DM_Network *network = (DM_Network*) dm->data; 1363a4e85ca8SHong Zhang PetscInt eStart,eEnd,vStart,vEnd,rstart,nrows,*rows,localSize; 1364840c2264SHong Zhang PetscInt cstart,ncols,j,e,v; 136524121865SAdrian Maldonado PetscBool ghost,ghost_vc,ghost2,isNest; 1366a4e85ca8SHong Zhang Mat Juser; 1367bfbc38dcSHong Zhang PetscSection sectionGlobal; 1368447d78afSSatish Balay PetscInt nedges,*vptr=NULL,vc,*rows_v; /* suppress maybe-uninitialized warning */ 1369a4e85ca8SHong Zhang const PetscInt *edges,*cone; 13705cf7da58SHong Zhang MPI_Comm comm; 137124121865SAdrian Maldonado MatType mtype; 13725cf7da58SHong Zhang Vec vd_nz,vo_nz; 13735cf7da58SHong Zhang PetscInt *dnnz,*onnz; 13745cf7da58SHong Zhang PetscScalar *vdnz,*vonz; 13751ad426b7SHong Zhang 13761ad426b7SHong Zhang PetscFunctionBegin; 137724121865SAdrian Maldonado mtype = dm->mattype; 137824121865SAdrian Maldonado ierr = PetscStrcmp(mtype, MATNEST, &isNest);CHKERRQ(ierr); 137924121865SAdrian Maldonado 138024121865SAdrian Maldonado if (isNest) { 13810731d606SHong Zhang /* ierr = DMCreateMatrix_Network_Nest(); */ 138224121865SAdrian Maldonado PetscInt eDof, vDof; 138324121865SAdrian Maldonado Mat j11, j12, j21, j22, bA[2][2]; 138424121865SAdrian Maldonado ISLocalToGlobalMapping eISMap, vISMap; 138524121865SAdrian Maldonado 138624121865SAdrian Maldonado ierr = PetscObjectGetComm((PetscObject)dm,&comm);CHKERRQ(ierr); 138724121865SAdrian Maldonado ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 138824121865SAdrian Maldonado ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 138924121865SAdrian Maldonado 139024121865SAdrian Maldonado ierr = PetscSectionGetConstrainedStorageSize(network->edge.GlobalDofSection,&eDof);CHKERRQ(ierr); 139124121865SAdrian Maldonado ierr = PetscSectionGetConstrainedStorageSize(network->vertex.GlobalDofSection,&vDof);CHKERRQ(ierr); 139224121865SAdrian Maldonado 139301ad2aeeSHong Zhang ierr = MatCreate(comm, &j11);CHKERRQ(ierr); 139424121865SAdrian Maldonado ierr = MatSetSizes(j11, eDof, eDof, PETSC_DETERMINE, PETSC_DETERMINE);CHKERRQ(ierr); 139524121865SAdrian Maldonado ierr = MatSetType(j11, MATMPIAIJ);CHKERRQ(ierr); 139624121865SAdrian Maldonado 139701ad2aeeSHong Zhang ierr = MatCreate(comm, &j12);CHKERRQ(ierr); 139824121865SAdrian Maldonado ierr = MatSetSizes(j12, eDof, vDof, PETSC_DETERMINE ,PETSC_DETERMINE);CHKERRQ(ierr); 139924121865SAdrian Maldonado ierr = MatSetType(j12, MATMPIAIJ);CHKERRQ(ierr); 140024121865SAdrian Maldonado 140101ad2aeeSHong Zhang ierr = MatCreate(comm, &j21);CHKERRQ(ierr); 140224121865SAdrian Maldonado ierr = MatSetSizes(j21, vDof, eDof, PETSC_DETERMINE, PETSC_DETERMINE);CHKERRQ(ierr); 140324121865SAdrian Maldonado ierr = MatSetType(j21, MATMPIAIJ);CHKERRQ(ierr); 140424121865SAdrian Maldonado 140501ad2aeeSHong Zhang ierr = MatCreate(comm, &j22);CHKERRQ(ierr); 140624121865SAdrian Maldonado ierr = MatSetSizes(j22, vDof, vDof, PETSC_DETERMINE, PETSC_DETERMINE);CHKERRQ(ierr); 140724121865SAdrian Maldonado ierr = MatSetType(j22, MATMPIAIJ);CHKERRQ(ierr); 140824121865SAdrian Maldonado 14093f6a6bdaSHong Zhang bA[0][0] = j11; 14103f6a6bdaSHong Zhang bA[0][1] = j12; 14113f6a6bdaSHong Zhang bA[1][0] = j21; 14123f6a6bdaSHong Zhang bA[1][1] = j22; 141324121865SAdrian Maldonado 141424121865SAdrian Maldonado ierr = CreateSubGlobalToLocalMapping_private(network->edge.GlobalDofSection,network->edge.DofSection,&eISMap);CHKERRQ(ierr); 141524121865SAdrian Maldonado ierr = CreateSubGlobalToLocalMapping_private(network->vertex.GlobalDofSection,network->vertex.DofSection,&vISMap);CHKERRQ(ierr); 141624121865SAdrian Maldonado 141724121865SAdrian Maldonado ierr = MatSetLocalToGlobalMapping(j11,eISMap,eISMap);CHKERRQ(ierr); 141824121865SAdrian Maldonado ierr = MatSetLocalToGlobalMapping(j12,eISMap,vISMap);CHKERRQ(ierr); 141924121865SAdrian Maldonado ierr = MatSetLocalToGlobalMapping(j21,vISMap,eISMap);CHKERRQ(ierr); 142024121865SAdrian Maldonado ierr = MatSetLocalToGlobalMapping(j22,vISMap,vISMap);CHKERRQ(ierr); 142124121865SAdrian Maldonado 142224121865SAdrian Maldonado ierr = MatSetUp(j11);CHKERRQ(ierr); 142324121865SAdrian Maldonado ierr = MatSetUp(j12);CHKERRQ(ierr); 142424121865SAdrian Maldonado ierr = MatSetUp(j21);CHKERRQ(ierr); 142524121865SAdrian Maldonado ierr = MatSetUp(j22);CHKERRQ(ierr); 142624121865SAdrian Maldonado 142701ad2aeeSHong Zhang ierr = MatCreateNest(comm,2,NULL,2,NULL,&bA[0][0],J);CHKERRQ(ierr); 142824121865SAdrian Maldonado ierr = MatSetUp(*J);CHKERRQ(ierr); 142924121865SAdrian Maldonado ierr = MatNestSetVecType(*J,VECNEST);CHKERRQ(ierr); 143024121865SAdrian Maldonado ierr = MatDestroy(&j11);CHKERRQ(ierr); 143124121865SAdrian Maldonado ierr = MatDestroy(&j12);CHKERRQ(ierr); 143224121865SAdrian Maldonado ierr = MatDestroy(&j21);CHKERRQ(ierr); 143324121865SAdrian Maldonado ierr = MatDestroy(&j22);CHKERRQ(ierr); 143424121865SAdrian Maldonado 143524121865SAdrian Maldonado ierr = MatAssemblyBegin(*J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 143624121865SAdrian Maldonado ierr = MatAssemblyEnd(*J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 143724121865SAdrian Maldonado 143824121865SAdrian Maldonado /* Free structures */ 143924121865SAdrian Maldonado ierr = ISLocalToGlobalMappingDestroy(&eISMap);CHKERRQ(ierr); 144024121865SAdrian Maldonado ierr = ISLocalToGlobalMappingDestroy(&vISMap);CHKERRQ(ierr); 144124121865SAdrian Maldonado 144224121865SAdrian Maldonado PetscFunctionReturn(0); 144324121865SAdrian Maldonado } else if (!network->userEdgeJacobian && !network->userVertexJacobian) { 1444a4e85ca8SHong Zhang /* user does not provide Jacobian blocks */ 1445bfbc38dcSHong Zhang ierr = DMCreateMatrix(network->plex,J);CHKERRQ(ierr); 1446bfbc38dcSHong Zhang ierr = MatSetDM(*J,dm);CHKERRQ(ierr); 14471ad426b7SHong Zhang PetscFunctionReturn(0); 14481ad426b7SHong Zhang } 14491ad426b7SHong Zhang 1450bfbc38dcSHong Zhang ierr = MatCreate(PetscObjectComm((PetscObject)dm),J);CHKERRQ(ierr); 14512a945128SHong Zhang ierr = DMGetDefaultGlobalSection(network->plex,§ionGlobal);CHKERRQ(ierr); 1452bfbc38dcSHong Zhang ierr = PetscSectionGetConstrainedStorageSize(sectionGlobal,&localSize);CHKERRQ(ierr); 1453bfbc38dcSHong Zhang ierr = MatSetSizes(*J,localSize,localSize,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr); 14542a945128SHong Zhang 14552a945128SHong Zhang ierr = MatSetType(*J,MATAIJ);CHKERRQ(ierr); 14562a945128SHong Zhang ierr = MatSetFromOptions(*J);CHKERRQ(ierr); 145789898e50SHong Zhang 145889898e50SHong Zhang /* (1) Set matrix preallocation */ 145989898e50SHong Zhang /*------------------------------*/ 1460840c2264SHong Zhang ierr = PetscObjectGetComm((PetscObject)dm,&comm);CHKERRQ(ierr); 1461840c2264SHong Zhang ierr = VecCreate(comm,&vd_nz);CHKERRQ(ierr); 1462840c2264SHong Zhang ierr = VecSetSizes(vd_nz,localSize,PETSC_DECIDE);CHKERRQ(ierr); 1463840c2264SHong Zhang ierr = VecSetFromOptions(vd_nz);CHKERRQ(ierr); 1464840c2264SHong Zhang ierr = VecSet(vd_nz,0.0);CHKERRQ(ierr); 1465840c2264SHong Zhang ierr = VecDuplicate(vd_nz,&vo_nz);CHKERRQ(ierr); 1466840c2264SHong Zhang 146789898e50SHong Zhang /* Set preallocation for edges */ 146889898e50SHong Zhang /*-----------------------------*/ 1469840c2264SHong Zhang ierr = DMNetworkGetEdgeRange(dm,&eStart,&eEnd);CHKERRQ(ierr); 1470840c2264SHong Zhang 1471bdcb62a2SHong Zhang ierr = PetscMalloc1(localSize,&rows);CHKERRQ(ierr); 1472840c2264SHong Zhang for (e=eStart; e<eEnd; e++) { 1473840c2264SHong Zhang /* Get row indices */ 1474840c2264SHong Zhang ierr = DMNetworkGetVariableGlobalOffset(dm,e,&rstart);CHKERRQ(ierr); 1475840c2264SHong Zhang ierr = DMNetworkGetNumVariables(dm,e,&nrows);CHKERRQ(ierr); 1476840c2264SHong Zhang if (nrows) { 1477840c2264SHong Zhang if (!network->Je) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_NULL,"Uer must provide Je"); 1478840c2264SHong Zhang for (j=0; j<nrows; j++) rows[j] = j + rstart; 1479840c2264SHong Zhang 14805cf7da58SHong Zhang /* Set preallocation for conntected vertices */ 1481d842c372SHong Zhang ierr = DMNetworkGetConnectedVertices(dm,e,&cone);CHKERRQ(ierr); 1482840c2264SHong Zhang for (v=0; v<2; v++) { 1483840c2264SHong Zhang ierr = DMNetworkGetNumVariables(dm,cone[v],&ncols);CHKERRQ(ierr); 1484840c2264SHong Zhang 1485840c2264SHong Zhang Juser = network->Je[3*e+1+v]; /* Jacobian(e,v) */ 1486840c2264SHong Zhang ierr = DMNetworkIsGhostVertex(dm,cone[v],&ghost);CHKERRQ(ierr); 14875cf7da58SHong Zhang ierr = MatSetPreallocationblock_private(Juser,nrows,rows,ncols,ghost,vd_nz,vo_nz);CHKERRQ(ierr); 1488840c2264SHong Zhang } 1489840c2264SHong Zhang 149089898e50SHong Zhang /* Set preallocation for edge self */ 1491840c2264SHong Zhang cstart = rstart; 1492840c2264SHong Zhang Juser = network->Je[3*e]; /* Jacobian(e,e) */ 14935cf7da58SHong Zhang ierr = MatSetPreallocationblock_private(Juser,nrows,rows,nrows,PETSC_FALSE,vd_nz,vo_nz);CHKERRQ(ierr); 1494840c2264SHong Zhang } 1495840c2264SHong Zhang } 1496840c2264SHong Zhang 149789898e50SHong Zhang /* Set preallocation for vertices */ 149889898e50SHong Zhang /*--------------------------------*/ 1499840c2264SHong Zhang ierr = DMNetworkGetVertexRange(dm,&vStart,&vEnd);CHKERRQ(ierr); 1500840c2264SHong Zhang if (vEnd - vStart) { 1501840c2264SHong Zhang if (!network->Jv) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_NULL,"Uer must provide Jv"); 1502840c2264SHong Zhang vptr = network->Jvptr; 1503840c2264SHong Zhang if (!vptr) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_NULL,"Uer must provide vptr"); 1504840c2264SHong Zhang } 1505840c2264SHong Zhang 1506840c2264SHong Zhang for (v=vStart; v<vEnd; v++) { 1507840c2264SHong Zhang /* Get row indices */ 1508840c2264SHong Zhang ierr = DMNetworkGetVariableGlobalOffset(dm,v,&rstart);CHKERRQ(ierr); 1509840c2264SHong Zhang ierr = DMNetworkGetNumVariables(dm,v,&nrows);CHKERRQ(ierr); 1510840c2264SHong Zhang if (!nrows) continue; 1511840c2264SHong Zhang 1512bdcb62a2SHong Zhang ierr = DMNetworkIsGhostVertex(dm,v,&ghost);CHKERRQ(ierr); 1513bdcb62a2SHong Zhang if (ghost) { 1514bdcb62a2SHong Zhang ierr = PetscMalloc1(nrows,&rows_v);CHKERRQ(ierr); 1515bdcb62a2SHong Zhang } else { 1516bdcb62a2SHong Zhang rows_v = rows; 1517bdcb62a2SHong Zhang } 1518bdcb62a2SHong Zhang 1519bdcb62a2SHong Zhang for (j=0; j<nrows; j++) rows_v[j] = j + rstart; 1520840c2264SHong Zhang 1521840c2264SHong Zhang /* Get supporting edges and connected vertices */ 1522840c2264SHong Zhang ierr = DMNetworkGetSupportingEdges(dm,v,&nedges,&edges);CHKERRQ(ierr); 1523840c2264SHong Zhang 1524840c2264SHong Zhang for (e=0; e<nedges; e++) { 1525840c2264SHong Zhang /* Supporting edges */ 1526840c2264SHong Zhang ierr = DMNetworkGetVariableGlobalOffset(dm,edges[e],&cstart);CHKERRQ(ierr); 1527840c2264SHong Zhang ierr = DMNetworkGetNumVariables(dm,edges[e],&ncols);CHKERRQ(ierr); 1528840c2264SHong Zhang 1529840c2264SHong Zhang Juser = network->Jv[vptr[v-vStart]+2*e+1]; /* Jacobian(v,e) */ 1530bdcb62a2SHong Zhang ierr = MatSetPreallocationblock_private(Juser,nrows,rows_v,ncols,ghost,vd_nz,vo_nz);CHKERRQ(ierr); 1531840c2264SHong Zhang 1532840c2264SHong Zhang /* Connected vertices */ 1533d842c372SHong Zhang ierr = DMNetworkGetConnectedVertices(dm,edges[e],&cone);CHKERRQ(ierr); 1534840c2264SHong Zhang vc = (v == cone[0]) ? cone[1]:cone[0]; 1535840c2264SHong Zhang ierr = DMNetworkIsGhostVertex(dm,vc,&ghost_vc);CHKERRQ(ierr); 1536840c2264SHong Zhang 1537840c2264SHong Zhang ierr = DMNetworkGetNumVariables(dm,vc,&ncols);CHKERRQ(ierr); 1538840c2264SHong Zhang 1539840c2264SHong Zhang Juser = network->Jv[vptr[v-vStart]+2*e+2]; /* Jacobian(v,vc) */ 1540e102a522SHong Zhang if (ghost_vc||ghost) { 1541e102a522SHong Zhang ghost2 = PETSC_TRUE; 1542e102a522SHong Zhang } else { 1543e102a522SHong Zhang ghost2 = PETSC_FALSE; 1544e102a522SHong Zhang } 1545e102a522SHong Zhang ierr = MatSetPreallocationblock_private(Juser,nrows,rows_v,ncols,ghost2,vd_nz,vo_nz);CHKERRQ(ierr); 1546840c2264SHong Zhang } 1547840c2264SHong Zhang 154889898e50SHong Zhang /* Set preallocation for vertex self */ 1549840c2264SHong Zhang ierr = DMNetworkIsGhostVertex(dm,v,&ghost);CHKERRQ(ierr); 1550840c2264SHong Zhang if (!ghost) { 1551840c2264SHong Zhang ierr = DMNetworkGetVariableGlobalOffset(dm,v,&cstart);CHKERRQ(ierr); 1552840c2264SHong Zhang Juser = network->Jv[vptr[v-vStart]]; /* Jacobian(v,v) */ 1553bdcb62a2SHong Zhang ierr = MatSetPreallocationblock_private(Juser,nrows,rows_v,nrows,PETSC_FALSE,vd_nz,vo_nz);CHKERRQ(ierr); 1554840c2264SHong Zhang } 1555bdcb62a2SHong Zhang if (ghost) { 1556bdcb62a2SHong Zhang ierr = PetscFree(rows_v);CHKERRQ(ierr); 1557bdcb62a2SHong Zhang } 1558840c2264SHong Zhang } 1559840c2264SHong Zhang 1560840c2264SHong Zhang ierr = VecAssemblyBegin(vd_nz);CHKERRQ(ierr); 1561840c2264SHong Zhang ierr = VecAssemblyBegin(vo_nz);CHKERRQ(ierr); 15625cf7da58SHong Zhang 15635cf7da58SHong Zhang ierr = PetscMalloc2(localSize,&dnnz,localSize,&onnz);CHKERRQ(ierr); 15645cf7da58SHong Zhang 15655cf7da58SHong Zhang ierr = VecAssemblyEnd(vd_nz);CHKERRQ(ierr); 1566840c2264SHong Zhang ierr = VecAssemblyEnd(vo_nz);CHKERRQ(ierr); 1567840c2264SHong Zhang 1568840c2264SHong Zhang ierr = VecGetArray(vd_nz,&vdnz);CHKERRQ(ierr); 1569840c2264SHong Zhang ierr = VecGetArray(vo_nz,&vonz);CHKERRQ(ierr); 1570840c2264SHong Zhang for (j=0; j<localSize; j++) { 1571e102a522SHong Zhang dnnz[j] = (PetscInt)PetscRealPart(vdnz[j]); 1572e102a522SHong Zhang onnz[j] = (PetscInt)PetscRealPart(vonz[j]); 1573840c2264SHong Zhang } 1574840c2264SHong Zhang ierr = VecRestoreArray(vd_nz,&vdnz);CHKERRQ(ierr); 1575840c2264SHong Zhang ierr = VecRestoreArray(vo_nz,&vonz);CHKERRQ(ierr); 1576840c2264SHong Zhang ierr = VecDestroy(&vd_nz);CHKERRQ(ierr); 1577840c2264SHong Zhang ierr = VecDestroy(&vo_nz);CHKERRQ(ierr); 1578840c2264SHong Zhang 15795cf7da58SHong Zhang ierr = MatSeqAIJSetPreallocation(*J,0,dnnz);CHKERRQ(ierr); 15805cf7da58SHong Zhang ierr = MatMPIAIJSetPreallocation(*J,0,dnnz,0,onnz);CHKERRQ(ierr); 15815cf7da58SHong Zhang ierr = MatSetOption(*J,MAT_NEW_NONZERO_ALLOCATION_ERR,PETSC_TRUE);CHKERRQ(ierr); 15825cf7da58SHong Zhang 15835cf7da58SHong Zhang ierr = PetscFree2(dnnz,onnz);CHKERRQ(ierr); 15845cf7da58SHong Zhang 158589898e50SHong Zhang /* (2) Set matrix entries for edges */ 158689898e50SHong Zhang /*----------------------------------*/ 15871ad426b7SHong Zhang for (e=eStart; e<eEnd; e++) { 1588bfbc38dcSHong Zhang /* Get row indices */ 15891ad426b7SHong Zhang ierr = DMNetworkGetVariableGlobalOffset(dm,e,&rstart);CHKERRQ(ierr); 159017df6e9eSHong Zhang ierr = DMNetworkGetNumVariables(dm,e,&nrows);CHKERRQ(ierr); 15914b976069SHong Zhang if (nrows) { 15924b976069SHong Zhang if (!network->Je) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_NULL,"Uer must provide Je"); 1593bdcb62a2SHong Zhang 159417df6e9eSHong Zhang for (j=0; j<nrows; j++) rows[j] = j + rstart; 15951ad426b7SHong Zhang 1596bfbc38dcSHong Zhang /* Set matrix entries for conntected vertices */ 1597d842c372SHong Zhang ierr = DMNetworkGetConnectedVertices(dm,e,&cone);CHKERRQ(ierr); 1598bfbc38dcSHong Zhang for (v=0; v<2; v++) { 1599bfbc38dcSHong Zhang ierr = DMNetworkGetVariableGlobalOffset(dm,cone[v],&cstart);CHKERRQ(ierr); 1600883e35e8SHong Zhang ierr = DMNetworkGetNumVariables(dm,cone[v],&ncols);CHKERRQ(ierr); 16013e97b6e8SHong Zhang 1602a4e85ca8SHong Zhang Juser = network->Je[3*e+1+v]; /* Jacobian(e,v) */ 1603a4e85ca8SHong Zhang ierr = MatSetblock_private(Juser,nrows,rows,ncols,cstart,J);CHKERRQ(ierr); 1604bfbc38dcSHong Zhang } 160517df6e9eSHong Zhang 1606bfbc38dcSHong Zhang /* Set matrix entries for edge self */ 16073e97b6e8SHong Zhang cstart = rstart; 1608a4e85ca8SHong Zhang Juser = network->Je[3*e]; /* Jacobian(e,e) */ 1609a4e85ca8SHong Zhang ierr = MatSetblock_private(Juser,nrows,rows,nrows,cstart,J);CHKERRQ(ierr); 16101ad426b7SHong Zhang } 16114b976069SHong Zhang } 16121ad426b7SHong Zhang 1613bfbc38dcSHong Zhang /* Set matrix entries for vertices */ 161483b2e829SHong Zhang /*---------------------------------*/ 16151ad426b7SHong Zhang for (v=vStart; v<vEnd; v++) { 1616bfbc38dcSHong Zhang /* Get row indices */ 1617596e729fSHong Zhang ierr = DMNetworkGetVariableGlobalOffset(dm,v,&rstart);CHKERRQ(ierr); 1618596e729fSHong Zhang ierr = DMNetworkGetNumVariables(dm,v,&nrows);CHKERRQ(ierr); 16194b976069SHong Zhang if (!nrows) continue; 1620596e729fSHong Zhang 1621bdcb62a2SHong Zhang ierr = DMNetworkIsGhostVertex(dm,v,&ghost);CHKERRQ(ierr); 1622bdcb62a2SHong Zhang if (ghost) { 1623bdcb62a2SHong Zhang ierr = PetscMalloc1(nrows,&rows_v);CHKERRQ(ierr); 1624bdcb62a2SHong Zhang } else { 1625bdcb62a2SHong Zhang rows_v = rows; 1626bdcb62a2SHong Zhang } 1627bdcb62a2SHong Zhang for (j=0; j<nrows; j++) rows_v[j] = j + rstart; 1628596e729fSHong Zhang 1629bfbc38dcSHong Zhang /* Get supporting edges and connected vertices */ 1630596e729fSHong Zhang ierr = DMNetworkGetSupportingEdges(dm,v,&nedges,&edges);CHKERRQ(ierr); 1631596e729fSHong Zhang 1632596e729fSHong Zhang for (e=0; e<nedges; e++) { 1633bfbc38dcSHong Zhang /* Supporting edges */ 1634596e729fSHong Zhang ierr = DMNetworkGetVariableGlobalOffset(dm,edges[e],&cstart);CHKERRQ(ierr); 1635596e729fSHong Zhang ierr = DMNetworkGetNumVariables(dm,edges[e],&ncols);CHKERRQ(ierr); 1636596e729fSHong Zhang 1637a4e85ca8SHong Zhang Juser = network->Jv[vptr[v-vStart]+2*e+1]; /* Jacobian(v,e) */ 1638bdcb62a2SHong Zhang ierr = MatSetblock_private(Juser,nrows,rows_v,ncols,cstart,J);CHKERRQ(ierr); 1639596e729fSHong Zhang 1640bfbc38dcSHong Zhang /* Connected vertices */ 1641d842c372SHong Zhang ierr = DMNetworkGetConnectedVertices(dm,edges[e],&cone);CHKERRQ(ierr); 16422a945128SHong Zhang vc = (v == cone[0]) ? cone[1]:cone[0]; 16432a945128SHong Zhang 164444aca652SHong Zhang ierr = DMNetworkGetVariableGlobalOffset(dm,vc,&cstart);CHKERRQ(ierr); 164544aca652SHong Zhang ierr = DMNetworkGetNumVariables(dm,vc,&ncols);CHKERRQ(ierr); 1646a4e85ca8SHong Zhang 1647a4e85ca8SHong Zhang Juser = network->Jv[vptr[v-vStart]+2*e+2]; /* Jacobian(v,vc) */ 1648bdcb62a2SHong Zhang ierr = MatSetblock_private(Juser,nrows,rows_v,ncols,cstart,J);CHKERRQ(ierr); 1649596e729fSHong Zhang } 1650596e729fSHong Zhang 1651bfbc38dcSHong Zhang /* Set matrix entries for vertex self */ 16521ad426b7SHong Zhang if (!ghost) { 1653596e729fSHong Zhang ierr = DMNetworkGetVariableGlobalOffset(dm,v,&cstart);CHKERRQ(ierr); 1654a4e85ca8SHong Zhang Juser = network->Jv[vptr[v-vStart]]; /* Jacobian(v,v) */ 1655bdcb62a2SHong Zhang ierr = MatSetblock_private(Juser,nrows,rows_v,nrows,cstart,J);CHKERRQ(ierr); 1656bdcb62a2SHong Zhang } 1657bdcb62a2SHong Zhang if (ghost) { 1658bdcb62a2SHong Zhang ierr = PetscFree(rows_v);CHKERRQ(ierr); 1659bdcb62a2SHong Zhang } 16601ad426b7SHong Zhang } 1661a4e85ca8SHong Zhang ierr = PetscFree(rows);CHKERRQ(ierr); 1662bdcb62a2SHong Zhang 16631ad426b7SHong Zhang ierr = MatAssemblyBegin(*J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 16641ad426b7SHong Zhang ierr = MatAssemblyEnd(*J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1665dd6f46cdSHong Zhang 16665f2c45f1SShri Abhyankar ierr = MatSetDM(*J,dm);CHKERRQ(ierr); 16675f2c45f1SShri Abhyankar PetscFunctionReturn(0); 16685f2c45f1SShri Abhyankar } 16695f2c45f1SShri Abhyankar 16705f2c45f1SShri Abhyankar PetscErrorCode DMDestroy_Network(DM dm) 16715f2c45f1SShri Abhyankar { 16725f2c45f1SShri Abhyankar PetscErrorCode ierr; 16735f2c45f1SShri Abhyankar DM_Network *network = (DM_Network*) dm->data; 1674*2727e31bSShri Abhyankar PetscInt j; 16755f2c45f1SShri Abhyankar 16765f2c45f1SShri Abhyankar PetscFunctionBegin; 16778415c774SShri Abhyankar if (--network->refct > 0) PetscFunctionReturn(0); 167883b2e829SHong Zhang if (network->Je) { 167983b2e829SHong Zhang ierr = PetscFree(network->Je);CHKERRQ(ierr); 168083b2e829SHong Zhang } 168183b2e829SHong Zhang if (network->Jv) { 1682883e35e8SHong Zhang ierr = PetscFree(network->Jvptr);CHKERRQ(ierr); 168383b2e829SHong Zhang ierr = PetscFree(network->Jv);CHKERRQ(ierr); 16841ad426b7SHong Zhang } 168513c2a604SAdrian Maldonado 168613c2a604SAdrian Maldonado ierr = ISLocalToGlobalMappingDestroy(&network->vertex.mapping);CHKERRQ(ierr); 168713c2a604SAdrian Maldonado ierr = PetscSectionDestroy(&network->vertex.DofSection);CHKERRQ(ierr); 168813c2a604SAdrian Maldonado ierr = PetscSectionDestroy(&network->vertex.GlobalDofSection);CHKERRQ(ierr); 168913c2a604SAdrian Maldonado if (network->vertex.sf) { 169013c2a604SAdrian Maldonado ierr = PetscSFDestroy(&network->vertex.sf);CHKERRQ(ierr); 169113c2a604SAdrian Maldonado } 169213c2a604SAdrian Maldonado /* edge */ 169313c2a604SAdrian Maldonado ierr = ISLocalToGlobalMappingDestroy(&network->edge.mapping);CHKERRQ(ierr); 169413c2a604SAdrian Maldonado ierr = PetscSectionDestroy(&network->edge.DofSection);CHKERRQ(ierr); 169513c2a604SAdrian Maldonado ierr = PetscSectionDestroy(&network->edge.GlobalDofSection);CHKERRQ(ierr); 169613c2a604SAdrian Maldonado if (network->edge.sf) { 169713c2a604SAdrian Maldonado ierr = PetscSFDestroy(&network->edge.sf);CHKERRQ(ierr); 169813c2a604SAdrian Maldonado } 16995f2c45f1SShri Abhyankar ierr = DMDestroy(&network->plex);CHKERRQ(ierr); 17005f2c45f1SShri Abhyankar network->edges = NULL; 17015f2c45f1SShri Abhyankar ierr = PetscSectionDestroy(&network->DataSection);CHKERRQ(ierr); 17025f2c45f1SShri Abhyankar ierr = PetscSectionDestroy(&network->DofSection);CHKERRQ(ierr); 170383b2e829SHong Zhang 1704*2727e31bSShri Abhyankar for(j=0; j < network->nsubnet; j++) { 1705*2727e31bSShri Abhyankar ierr = PetscFree(network->subnet[j].edges);CHKERRQ(ierr); 1706*2727e31bSShri Abhyankar ierr = PetscFree(network->subnet[j].vertices);CHKERRQ(ierr); 1707*2727e31bSShri Abhyankar } 1708e2aaf10cSShri Abhyankar ierr = PetscFree(network->subnet);CHKERRQ(ierr); 17095f2c45f1SShri Abhyankar ierr = PetscFree(network->componentdataarray);CHKERRQ(ierr); 17105f2c45f1SShri Abhyankar ierr = PetscFree(network->cvalue);CHKERRQ(ierr); 17115f2c45f1SShri Abhyankar ierr = PetscFree(network->header);CHKERRQ(ierr); 17125f2c45f1SShri Abhyankar ierr = PetscFree(network);CHKERRQ(ierr); 17135f2c45f1SShri Abhyankar PetscFunctionReturn(0); 17145f2c45f1SShri Abhyankar } 17155f2c45f1SShri Abhyankar 17165f2c45f1SShri Abhyankar PetscErrorCode DMView_Network(DM dm, PetscViewer viewer) 17175f2c45f1SShri Abhyankar { 17185f2c45f1SShri Abhyankar PetscErrorCode ierr; 17195f2c45f1SShri Abhyankar DM_Network *network = (DM_Network*) dm->data; 17205f2c45f1SShri Abhyankar 17215f2c45f1SShri Abhyankar PetscFunctionBegin; 17225f2c45f1SShri Abhyankar ierr = DMView(network->plex,viewer);CHKERRQ(ierr); 17235f2c45f1SShri Abhyankar PetscFunctionReturn(0); 17245f2c45f1SShri Abhyankar } 17255f2c45f1SShri Abhyankar 17265f2c45f1SShri Abhyankar PetscErrorCode DMGlobalToLocalBegin_Network(DM dm, Vec g, InsertMode mode, Vec l) 17275f2c45f1SShri Abhyankar { 17285f2c45f1SShri Abhyankar PetscErrorCode ierr; 17295f2c45f1SShri Abhyankar DM_Network *network = (DM_Network*) dm->data; 17305f2c45f1SShri Abhyankar 17315f2c45f1SShri Abhyankar PetscFunctionBegin; 17325f2c45f1SShri Abhyankar ierr = DMGlobalToLocalBegin(network->plex,g,mode,l);CHKERRQ(ierr); 17335f2c45f1SShri Abhyankar PetscFunctionReturn(0); 17345f2c45f1SShri Abhyankar } 17355f2c45f1SShri Abhyankar 17365f2c45f1SShri Abhyankar PetscErrorCode DMGlobalToLocalEnd_Network(DM dm, Vec g, InsertMode mode, Vec l) 17375f2c45f1SShri Abhyankar { 17385f2c45f1SShri Abhyankar PetscErrorCode ierr; 17395f2c45f1SShri Abhyankar DM_Network *network = (DM_Network*) dm->data; 17405f2c45f1SShri Abhyankar 17415f2c45f1SShri Abhyankar PetscFunctionBegin; 17425f2c45f1SShri Abhyankar ierr = DMGlobalToLocalEnd(network->plex,g,mode,l);CHKERRQ(ierr); 17435f2c45f1SShri Abhyankar PetscFunctionReturn(0); 17445f2c45f1SShri Abhyankar } 17455f2c45f1SShri Abhyankar 17465f2c45f1SShri Abhyankar PetscErrorCode DMLocalToGlobalBegin_Network(DM dm, Vec l, InsertMode mode, Vec g) 17475f2c45f1SShri Abhyankar { 17485f2c45f1SShri Abhyankar PetscErrorCode ierr; 17495f2c45f1SShri Abhyankar DM_Network *network = (DM_Network*) dm->data; 17505f2c45f1SShri Abhyankar 17515f2c45f1SShri Abhyankar PetscFunctionBegin; 17525f2c45f1SShri Abhyankar ierr = DMLocalToGlobalBegin(network->plex,l,mode,g);CHKERRQ(ierr); 17535f2c45f1SShri Abhyankar PetscFunctionReturn(0); 17545f2c45f1SShri Abhyankar } 17555f2c45f1SShri Abhyankar 17565f2c45f1SShri Abhyankar PetscErrorCode DMLocalToGlobalEnd_Network(DM dm, Vec l, InsertMode mode, Vec g) 17575f2c45f1SShri Abhyankar { 17585f2c45f1SShri Abhyankar PetscErrorCode ierr; 17595f2c45f1SShri Abhyankar DM_Network *network = (DM_Network*) dm->data; 17605f2c45f1SShri Abhyankar 17615f2c45f1SShri Abhyankar PetscFunctionBegin; 17625f2c45f1SShri Abhyankar ierr = DMLocalToGlobalEnd(network->plex,l,mode,g);CHKERRQ(ierr); 17635f2c45f1SShri Abhyankar PetscFunctionReturn(0); 17645f2c45f1SShri Abhyankar } 1765