1*5f2c45f1SShri Abhyankar #include <petsc-private/dmnetworkimpl.h> /*I "petscdmnetwork.h" I*/ 2*5f2c45f1SShri Abhyankar #include <petscdmplex.h> 3*5f2c45f1SShri Abhyankar #include <petscsf.h> 4*5f2c45f1SShri Abhyankar 5*5f2c45f1SShri Abhyankar #undef __FUNCT__ 6*5f2c45f1SShri Abhyankar #define __FUNCT__ "DMNetworkSetSizes" 7*5f2c45f1SShri Abhyankar /*@ 8*5f2c45f1SShri Abhyankar DMNetworkSetSizes - Sets the local and global vertices and edges. 9*5f2c45f1SShri Abhyankar 10*5f2c45f1SShri Abhyankar Collective on DM 11*5f2c45f1SShri Abhyankar 12*5f2c45f1SShri Abhyankar Input Parameters: 13*5f2c45f1SShri Abhyankar + dm - the dm object 14*5f2c45f1SShri Abhyankar . nV - number of local vertices 15*5f2c45f1SShri Abhyankar . nE - number of local edges 16*5f2c45f1SShri Abhyankar . NV - number of global vertices (or PETSC_DETERMINE) 17*5f2c45f1SShri Abhyankar - NE - number of global edges (or PETSC_DETERMINE) 18*5f2c45f1SShri Abhyankar 19*5f2c45f1SShri Abhyankar Notes 20*5f2c45f1SShri Abhyankar If one processor calls this with NV (NE) of PETSC_DECIDE then all processors must, otherwise the prgram will hang. 21*5f2c45f1SShri Abhyankar 22*5f2c45f1SShri Abhyankar You cannot change the sizes once they have been set 23*5f2c45f1SShri Abhyankar 24*5f2c45f1SShri Abhyankar Level: intermediate 25*5f2c45f1SShri Abhyankar 26*5f2c45f1SShri Abhyankar .seealso: DMNetworkCreate 27*5f2c45f1SShri Abhyankar @*/ 28*5f2c45f1SShri Abhyankar PetscErrorCode DMNetworkSetSizes(DM dm, PetscInt nV, PetscInt nE, PetscInt NV, PetscInt NE) 29*5f2c45f1SShri Abhyankar { 30*5f2c45f1SShri Abhyankar PetscErrorCode ierr; 31*5f2c45f1SShri Abhyankar DM_Network *network = (DM_Network*) dm->data; 32*5f2c45f1SShri Abhyankar PetscInt a[2],b[2]; 33*5f2c45f1SShri Abhyankar 34*5f2c45f1SShri Abhyankar PetscFunctionBegin; 35*5f2c45f1SShri Abhyankar PetscValidHeaderSpecific(dm,DM_CLASSID,1); 36*5f2c45f1SShri Abhyankar if (NV > 0) PetscValidLogicalCollectiveInt(dm,NV,4); 37*5f2c45f1SShri Abhyankar if (NE > 0) PetscValidLogicalCollectiveInt(dm,NE,5); 38*5f2c45f1SShri Abhyankar if (NV > 0 && nV > NV) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Local vertex size %D cannot be larger than global vertex size %D",nV,NV); 39*5f2c45f1SShri Abhyankar if (NE > 0 && nE > NE) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Local edge size %D cannot be larger than global edge size %D",nE,NE); 40*5f2c45f1SShri Abhyankar if ((network->nNodes >= 0 || network->NNodes >= 0) && (network->nNodes != nV || network->NNodes != NV)) SETERRQ4(PETSC_COMM_SELF,PETSC_ERR_SUP,"Cannot change/reset vertex sizes to %D local %D global after previously setting them to %D local %D global",nV,NV,network->nNodes,network->NNodes); 41*5f2c45f1SShri Abhyankar if ((network->nEdges >= 0 || network->NEdges >= 0) && (network->nEdges != nE || network->NEdges != NE)) SETERRQ4(PETSC_COMM_SELF,PETSC_ERR_SUP,"Cannot change/reset edge sizes to %D local %D global after previously setting them to %D local %D global",nE,NE,network->nEdges,network->NEdges); 42*5f2c45f1SShri Abhyankar if (NE < 0 || NV < 0) { 43*5f2c45f1SShri Abhyankar a[0] = nV; a[1] = nE; 44*5f2c45f1SShri Abhyankar ierr = MPI_Allreduce(a,b,2,MPIU_INT,MPI_SUM,PetscObjectComm((PetscObject)dm));CHKERRQ(ierr); 45*5f2c45f1SShri Abhyankar NV = b[0]; NE = b[1]; 46*5f2c45f1SShri Abhyankar } 47*5f2c45f1SShri Abhyankar network->nNodes = nV; 48*5f2c45f1SShri Abhyankar network->NNodes = NV; 49*5f2c45f1SShri Abhyankar network->nEdges = nE; 50*5f2c45f1SShri Abhyankar network->NEdges = NE; 51*5f2c45f1SShri Abhyankar PetscFunctionReturn(0); 52*5f2c45f1SShri Abhyankar } 53*5f2c45f1SShri Abhyankar 54*5f2c45f1SShri Abhyankar #undef __FUNCT__ 55*5f2c45f1SShri Abhyankar #define __FUNCT__ "DMNetworkSetEdgeList" 56*5f2c45f1SShri Abhyankar /*@ 57*5f2c45f1SShri Abhyankar DMNetworkSetEdgeList - Sets the list of local edges (vertex connectivity) for the network 58*5f2c45f1SShri Abhyankar 59*5f2c45f1SShri Abhyankar Logically collective on DM 60*5f2c45f1SShri Abhyankar 61*5f2c45f1SShri Abhyankar Input Parameters: 62*5f2c45f1SShri Abhyankar . edges - list of edges 63*5f2c45f1SShri Abhyankar 64*5f2c45f1SShri Abhyankar Notes: 65*5f2c45f1SShri Abhyankar There is no copy involved in this operation, only the pointer is referenced. The edgelist should 66*5f2c45f1SShri Abhyankar not be destroyed before the call to DMNetworkLayoutSetUp 67*5f2c45f1SShri Abhyankar 68*5f2c45f1SShri Abhyankar Level: intermediate 69*5f2c45f1SShri Abhyankar 70*5f2c45f1SShri Abhyankar .seealso: DMNetworkCreate, DMNetworkSetSizes 71*5f2c45f1SShri Abhyankar @*/ 72*5f2c45f1SShri Abhyankar PetscErrorCode DMNetworkSetEdgeList(DM dm, int edgelist[]) 73*5f2c45f1SShri Abhyankar { 74*5f2c45f1SShri Abhyankar DM_Network *network = (DM_Network*) dm->data; 75*5f2c45f1SShri Abhyankar 76*5f2c45f1SShri Abhyankar PetscFunctionBegin; 77*5f2c45f1SShri Abhyankar network->edges = edgelist; 78*5f2c45f1SShri Abhyankar PetscFunctionReturn(0); 79*5f2c45f1SShri Abhyankar } 80*5f2c45f1SShri Abhyankar 81*5f2c45f1SShri Abhyankar #undef __FUNCT__ 82*5f2c45f1SShri Abhyankar #define __FUNCT__ "DMNetworkLayoutSetUp" 83*5f2c45f1SShri Abhyankar /*@ 84*5f2c45f1SShri Abhyankar DMNetworkLayoutSetUp - Sets up the bare layout (graph) for the network 85*5f2c45f1SShri Abhyankar 86*5f2c45f1SShri Abhyankar Collective on DM 87*5f2c45f1SShri Abhyankar 88*5f2c45f1SShri Abhyankar Input Parameters 89*5f2c45f1SShri Abhyankar . DM - the dmnetwork object 90*5f2c45f1SShri Abhyankar 91*5f2c45f1SShri Abhyankar Notes: 92*5f2c45f1SShri Abhyankar This routine should be called after the network sizes and edgelists have been provided. It creates 93*5f2c45f1SShri Abhyankar the bare layout of the network and sets up the network to begin insertion of components. 94*5f2c45f1SShri Abhyankar 95*5f2c45f1SShri Abhyankar All the components should be registered before calling this routine. 96*5f2c45f1SShri Abhyankar 97*5f2c45f1SShri Abhyankar Level: intermediate 98*5f2c45f1SShri Abhyankar 99*5f2c45f1SShri Abhyankar .seealso: DMNetworkSetSizes, DMNetworkSetEdgeList 100*5f2c45f1SShri Abhyankar @*/ 101*5f2c45f1SShri Abhyankar PetscErrorCode DMNetworkLayoutSetUp(DM dm) 102*5f2c45f1SShri Abhyankar { 103*5f2c45f1SShri Abhyankar PetscErrorCode ierr; 104*5f2c45f1SShri Abhyankar DM_Network *network = (DM_Network*) dm->data; 105*5f2c45f1SShri Abhyankar PetscInt dim = 1; /* One dimensional network */ 106*5f2c45f1SShri Abhyankar PetscInt numCorners=2; 107*5f2c45f1SShri Abhyankar PetscInt spacedim=2; 108*5f2c45f1SShri Abhyankar double *vertexcoords=NULL; 109*5f2c45f1SShri Abhyankar PetscInt i; 110*5f2c45f1SShri Abhyankar PetscInt ndata; 111*5f2c45f1SShri Abhyankar 112*5f2c45f1SShri Abhyankar PetscFunctionBegin; 113*5f2c45f1SShri Abhyankar if (network->nNodes) { 114*5f2c45f1SShri Abhyankar ierr = PetscMalloc(numCorners*network->nNodes*sizeof(PetscInt),&vertexcoords);CHKERRQ(ierr); 115*5f2c45f1SShri Abhyankar } 116*5f2c45f1SShri Abhyankar ierr = DMPlexCreateFromCellList(PetscObjectComm((PetscObject)dm),dim,network->nEdges,network->nNodes,numCorners,PETSC_FALSE,network->edges,spacedim,vertexcoords,&network->plex);CHKERRQ(ierr); 117*5f2c45f1SShri Abhyankar if (network->nNodes) { 118*5f2c45f1SShri Abhyankar ierr = PetscFree(vertexcoords);CHKERRQ(ierr); 119*5f2c45f1SShri Abhyankar } 120*5f2c45f1SShri Abhyankar ierr = DMPlexGetChart(network->plex,&network->pStart,&network->pEnd);CHKERRQ(ierr); 121*5f2c45f1SShri Abhyankar ierr = DMPlexGetHeightStratum(network->plex,0,&network->eStart,&network->eEnd);CHKERRQ(ierr); 122*5f2c45f1SShri Abhyankar ierr = DMPlexGetHeightStratum(network->plex,1,&network->vStart,&network->vEnd);CHKERRQ(ierr); 123*5f2c45f1SShri Abhyankar 124*5f2c45f1SShri Abhyankar ierr = PetscSectionCreate(PetscObjectComm((PetscObject)dm),&network->DataSection);CHKERRQ(ierr); 125*5f2c45f1SShri Abhyankar ierr = PetscSectionCreate(PetscObjectComm((PetscObject)dm),&network->DofSection);CHKERRQ(ierr); 126*5f2c45f1SShri Abhyankar ierr = PetscSectionSetChart(network->DataSection,network->pStart,network->pEnd);CHKERRQ(ierr); 127*5f2c45f1SShri Abhyankar ierr = PetscSectionSetChart(network->DofSection,network->pStart,network->pEnd);CHKERRQ(ierr); 128*5f2c45f1SShri Abhyankar 129*5f2c45f1SShri Abhyankar network->dataheadersize = sizeof(struct _p_DMNetworkComponentHeader)/sizeof(DMNetworkComponentGenericDataType); 130*5f2c45f1SShri Abhyankar ierr = PetscMalloc((network->pEnd-network->pStart)*sizeof(struct _p_DMNetworkComponentHeader),&network->header);CHKERRQ(ierr); 131*5f2c45f1SShri Abhyankar for (i = network->pStart; i < network->pEnd; i++) { 132*5f2c45f1SShri Abhyankar network->header[i].ndata = 0; 133*5f2c45f1SShri Abhyankar ndata = network->header[i].ndata; 134*5f2c45f1SShri Abhyankar ierr = PetscSectionAddDof(network->DataSection,i,network->dataheadersize);CHKERRQ(ierr); 135*5f2c45f1SShri Abhyankar network->header[i].offset[ndata] = 0; 136*5f2c45f1SShri Abhyankar } 137*5f2c45f1SShri Abhyankar ierr = PetscMalloc((network->pEnd-network->pStart)*sizeof(struct _p_DMNetworkComponentValue),&network->cvalue);CHKERRQ(ierr); 138*5f2c45f1SShri Abhyankar PetscFunctionReturn(0); 139*5f2c45f1SShri Abhyankar } 140*5f2c45f1SShri Abhyankar 141*5f2c45f1SShri Abhyankar #undef __FUNCT__ 142*5f2c45f1SShri Abhyankar #define __FUNCT__ "DMNetworkRegisterComponent" 143*5f2c45f1SShri Abhyankar /*@ 144*5f2c45f1SShri Abhyankar DMNetworkRegisterComponent - Registers the network component 145*5f2c45f1SShri Abhyankar 146*5f2c45f1SShri Abhyankar Logically collective on DM 147*5f2c45f1SShri Abhyankar 148*5f2c45f1SShri Abhyankar Input Parameters 149*5f2c45f1SShri Abhyankar + dm - the network object 150*5f2c45f1SShri Abhyankar . name - the component name 151*5f2c45f1SShri Abhyankar - size - the storage size in bytes for this component data 152*5f2c45f1SShri Abhyankar 153*5f2c45f1SShri Abhyankar Output Parameters 154*5f2c45f1SShri Abhyankar . key - an integer key that defines the component 155*5f2c45f1SShri Abhyankar 156*5f2c45f1SShri Abhyankar Notes 157*5f2c45f1SShri Abhyankar This routine should be called by all processors before calling DMNetworkLayoutSetup(). 158*5f2c45f1SShri Abhyankar 159*5f2c45f1SShri Abhyankar Level: intermediate 160*5f2c45f1SShri Abhyankar 161*5f2c45f1SShri Abhyankar .seealso: DMNetworkLayoutSetUp, DMNetworkCreate 162*5f2c45f1SShri Abhyankar @*/ 163*5f2c45f1SShri Abhyankar PetscErrorCode DMNetworkRegisterComponent(DM dm,const char *name,PetscInt size,PetscInt *key) 164*5f2c45f1SShri Abhyankar { 165*5f2c45f1SShri Abhyankar PetscErrorCode ierr; 166*5f2c45f1SShri Abhyankar DM_Network *network = (DM_Network*) dm->data; 167*5f2c45f1SShri Abhyankar DMNetworkComponent *component=&network->component[network->ncomponent]; 168*5f2c45f1SShri Abhyankar PetscBool flg=PETSC_FALSE; 169*5f2c45f1SShri Abhyankar PetscInt i; 170*5f2c45f1SShri Abhyankar 171*5f2c45f1SShri Abhyankar PetscFunctionBegin; 172*5f2c45f1SShri Abhyankar 173*5f2c45f1SShri Abhyankar for (i=0; i < network->ncomponent; i++) { 174*5f2c45f1SShri Abhyankar ierr = PetscStrcmp(component->name,name,&flg);CHKERRQ(ierr); 175*5f2c45f1SShri Abhyankar if (flg) { 176*5f2c45f1SShri Abhyankar *key = i; 177*5f2c45f1SShri Abhyankar PetscFunctionReturn(0); 178*5f2c45f1SShri Abhyankar } 179*5f2c45f1SShri Abhyankar } 180*5f2c45f1SShri Abhyankar 181*5f2c45f1SShri Abhyankar ierr = PetscStrcpy(component->name,name);CHKERRQ(ierr); 182*5f2c45f1SShri Abhyankar component->size = size/sizeof(DMNetworkComponentGenericDataType); 183*5f2c45f1SShri Abhyankar *key = network->ncomponent; 184*5f2c45f1SShri Abhyankar network->ncomponent++; 185*5f2c45f1SShri Abhyankar PetscFunctionReturn(0); 186*5f2c45f1SShri Abhyankar } 187*5f2c45f1SShri Abhyankar 188*5f2c45f1SShri Abhyankar #undef __FUNCT__ 189*5f2c45f1SShri Abhyankar #define __FUNCT__ "DMNetworkGetVertexRange" 190*5f2c45f1SShri Abhyankar /*@ 191*5f2c45f1SShri Abhyankar DMNetworkGetVertexRange - Get the bounds [start, end) for the vertices. 192*5f2c45f1SShri Abhyankar 193*5f2c45f1SShri Abhyankar Not Collective 194*5f2c45f1SShri Abhyankar 195*5f2c45f1SShri Abhyankar Input Parameters: 196*5f2c45f1SShri Abhyankar + dm - The DMNetwork object 197*5f2c45f1SShri Abhyankar 198*5f2c45f1SShri Abhyankar Output Paramters: 199*5f2c45f1SShri Abhyankar + vStart - The first vertex point 200*5f2c45f1SShri Abhyankar - vEnd - One beyond the last vertex point 201*5f2c45f1SShri Abhyankar 202*5f2c45f1SShri Abhyankar Level: intermediate 203*5f2c45f1SShri Abhyankar 204*5f2c45f1SShri Abhyankar .seealso: DMNetworkGetEdgeRange 205*5f2c45f1SShri Abhyankar @*/ 206*5f2c45f1SShri Abhyankar PetscErrorCode DMNetworkGetVertexRange(DM dm,PetscInt *vStart,PetscInt *vEnd) 207*5f2c45f1SShri Abhyankar { 208*5f2c45f1SShri Abhyankar DM_Network *network = (DM_Network*)dm->data; 209*5f2c45f1SShri Abhyankar 210*5f2c45f1SShri Abhyankar PetscFunctionBegin; 211*5f2c45f1SShri Abhyankar if (vStart) *vStart = network->vStart; 212*5f2c45f1SShri Abhyankar if (vEnd) *vEnd = network->vEnd; 213*5f2c45f1SShri Abhyankar PetscFunctionReturn(0); 214*5f2c45f1SShri Abhyankar } 215*5f2c45f1SShri Abhyankar 216*5f2c45f1SShri Abhyankar #undef __FUNCT__ 217*5f2c45f1SShri Abhyankar #define __FUNCT__ "DMNetworkGetEdgeRange" 218*5f2c45f1SShri Abhyankar /*@ 219*5f2c45f1SShri Abhyankar DMNetworkGetEdgeRange - Get the bounds [start, end) for the edges. 220*5f2c45f1SShri Abhyankar 221*5f2c45f1SShri Abhyankar Not Collective 222*5f2c45f1SShri Abhyankar 223*5f2c45f1SShri Abhyankar Input Parameters: 224*5f2c45f1SShri Abhyankar + dm - The DMNetwork object 225*5f2c45f1SShri Abhyankar 226*5f2c45f1SShri Abhyankar Output Paramters: 227*5f2c45f1SShri Abhyankar + eStart - The first edge point 228*5f2c45f1SShri Abhyankar - eEnd - One beyond the last edge point 229*5f2c45f1SShri Abhyankar 230*5f2c45f1SShri Abhyankar Level: intermediate 231*5f2c45f1SShri Abhyankar 232*5f2c45f1SShri Abhyankar .seealso: DMNetworkGetVertexRange 233*5f2c45f1SShri Abhyankar @*/ 234*5f2c45f1SShri Abhyankar PetscErrorCode DMNetworkGetEdgeRange(DM dm,PetscInt *eStart,PetscInt *eEnd) 235*5f2c45f1SShri Abhyankar { 236*5f2c45f1SShri Abhyankar DM_Network *network = (DM_Network*)dm->data; 237*5f2c45f1SShri Abhyankar 238*5f2c45f1SShri Abhyankar PetscFunctionBegin; 239*5f2c45f1SShri Abhyankar if (eStart) *eStart = network->eStart; 240*5f2c45f1SShri Abhyankar if (eEnd) *eEnd = network->eEnd; 241*5f2c45f1SShri Abhyankar PetscFunctionReturn(0); 242*5f2c45f1SShri Abhyankar } 243*5f2c45f1SShri Abhyankar 244*5f2c45f1SShri Abhyankar #undef __FUNCT__ 245*5f2c45f1SShri Abhyankar #define __FUNCT__ "DMNetworkAddComponent" 246*5f2c45f1SShri Abhyankar /*@ 247*5f2c45f1SShri Abhyankar DMNetworkAddComponent - Adds a network component at the given point (vertex/edge) 248*5f2c45f1SShri Abhyankar 249*5f2c45f1SShri Abhyankar Not Collective 250*5f2c45f1SShri Abhyankar 251*5f2c45f1SShri Abhyankar Input Parameters: 252*5f2c45f1SShri Abhyankar + dm - The DMNetwork object 253*5f2c45f1SShri Abhyankar . p - vertex/edge point 254*5f2c45f1SShri Abhyankar . componentkey - component key returned while registering the component 255*5f2c45f1SShri Abhyankar - compvalue - pointer to the data structure for the component 256*5f2c45f1SShri Abhyankar 257*5f2c45f1SShri Abhyankar Level: intermediate 258*5f2c45f1SShri Abhyankar 259*5f2c45f1SShri Abhyankar .seealso: DMNetworkGetVertexRange, DMNetworkGetEdgeRange, DMNetworkRegisterComponent 260*5f2c45f1SShri Abhyankar @*/ 261*5f2c45f1SShri Abhyankar PetscErrorCode DMNetworkAddComponent(DM dm, PetscInt p,PetscInt componentkey,void* compvalue) 262*5f2c45f1SShri Abhyankar { 263*5f2c45f1SShri Abhyankar DM_Network *network = (DM_Network*)dm->data; 264*5f2c45f1SShri Abhyankar DMNetworkComponent component=network->component[componentkey]; 265*5f2c45f1SShri Abhyankar DMNetworkComponentHeader header=&network->header[p]; 266*5f2c45f1SShri Abhyankar DMNetworkComponentValue cvalue=&network->cvalue[p]; 267*5f2c45f1SShri Abhyankar PetscErrorCode ierr; 268*5f2c45f1SShri Abhyankar 269*5f2c45f1SShri Abhyankar PetscFunctionBegin; 270*5f2c45f1SShri Abhyankar header->size[header->ndata] = component.size; 271*5f2c45f1SShri Abhyankar ierr = PetscSectionAddDof(network->DataSection,p,component.size);CHKERRQ(ierr); 272*5f2c45f1SShri Abhyankar header->key[header->ndata] = componentkey; 273*5f2c45f1SShri Abhyankar if (header->ndata != 0) header->offset[header->ndata] = header->offset[header->ndata-1] + header->size[header->ndata-1]; 274*5f2c45f1SShri Abhyankar 275*5f2c45f1SShri Abhyankar cvalue->data[header->ndata] = (void*)compvalue; 276*5f2c45f1SShri Abhyankar header->ndata++; 277*5f2c45f1SShri Abhyankar PetscFunctionReturn(0); 278*5f2c45f1SShri Abhyankar } 279*5f2c45f1SShri Abhyankar 280*5f2c45f1SShri Abhyankar #undef __FUNCT__ 281*5f2c45f1SShri Abhyankar #define __FUNCT__ "DMNetworkGetNumComponents" 282*5f2c45f1SShri Abhyankar /*@ 283*5f2c45f1SShri Abhyankar DMNetworkGetNumComponents - Get the number of components at a vertex/edge 284*5f2c45f1SShri Abhyankar 285*5f2c45f1SShri Abhyankar Not Collective 286*5f2c45f1SShri Abhyankar 287*5f2c45f1SShri Abhyankar Input Parameters: 288*5f2c45f1SShri Abhyankar + dm - The DMNetwork object 289*5f2c45f1SShri Abhyankar . p - vertex/edge point 290*5f2c45f1SShri Abhyankar 291*5f2c45f1SShri Abhyankar Output Parameters: 292*5f2c45f1SShri Abhyankar . numcomponents - Number of components at the vertex/edge 293*5f2c45f1SShri Abhyankar 294*5f2c45f1SShri Abhyankar Level: intermediate 295*5f2c45f1SShri Abhyankar 296*5f2c45f1SShri Abhyankar .seealso: DMNetworkRegisterComponent, DMNetworkAddComponent 297*5f2c45f1SShri Abhyankar @*/ 298*5f2c45f1SShri Abhyankar PetscErrorCode DMNetworkGetNumComponents(DM dm,PetscInt p,PetscInt *numcomponents) 299*5f2c45f1SShri Abhyankar { 300*5f2c45f1SShri Abhyankar PetscErrorCode ierr; 301*5f2c45f1SShri Abhyankar PetscInt offset; 302*5f2c45f1SShri Abhyankar DM_Network *network = (DM_Network*)dm->data; 303*5f2c45f1SShri Abhyankar 304*5f2c45f1SShri Abhyankar PetscFunctionBegin; 305*5f2c45f1SShri Abhyankar ierr = PetscSectionGetOffset(network->DataSection,p,&offset);CHKERRQ(ierr); 306*5f2c45f1SShri Abhyankar *numcomponents = ((DMNetworkComponentHeader)(network->componentdataarray+offset))->ndata; 307*5f2c45f1SShri Abhyankar PetscFunctionReturn(0); 308*5f2c45f1SShri Abhyankar } 309*5f2c45f1SShri Abhyankar 310*5f2c45f1SShri Abhyankar #undef __FUNCT__ 311*5f2c45f1SShri Abhyankar #define __FUNCT__ "DMNetworkGetComponentTypeOffset" 312*5f2c45f1SShri Abhyankar /*@ 313*5f2c45f1SShri Abhyankar DMNetworkGetComponentTypeOffset - Gets the type along with the offset for indexing the 314*5f2c45f1SShri Abhyankar component value from the component data array 315*5f2c45f1SShri Abhyankar 316*5f2c45f1SShri Abhyankar Not Collective 317*5f2c45f1SShri Abhyankar 318*5f2c45f1SShri Abhyankar Input Parameters: 319*5f2c45f1SShri Abhyankar + dm - The DMNetwork object 320*5f2c45f1SShri Abhyankar . p - vertex/edge point 321*5f2c45f1SShri Abhyankar - compnum - component number 322*5f2c45f1SShri Abhyankar 323*5f2c45f1SShri Abhyankar Output Parameters: 324*5f2c45f1SShri Abhyankar + compkey - the key obtained when registering the component 325*5f2c45f1SShri Abhyankar - offset - offset into the component data array associated with the vertex/edge point 326*5f2c45f1SShri Abhyankar 327*5f2c45f1SShri Abhyankar Notes: 328*5f2c45f1SShri Abhyankar Typical usage: 329*5f2c45f1SShri Abhyankar 330*5f2c45f1SShri Abhyankar DMNetworkGetComponentDataArray(dm, &arr); 331*5f2c45f1SShri Abhyankar DMNetworkGetVertex/EdgeRange(dm,&Start,&End); 332*5f2c45f1SShri Abhyankar Loop over vertices or edges 333*5f2c45f1SShri Abhyankar DMNetworkGetNumComponents(dm,v,&numcomps); 334*5f2c45f1SShri Abhyankar Loop over numcomps 335*5f2c45f1SShri Abhyankar DMNetworkGetComponentTypeOffset(dm,v,compnum,&key,&offset); 336*5f2c45f1SShri Abhyankar compdata = (UserCompDataType)(arr+offset); 337*5f2c45f1SShri Abhyankar 338*5f2c45f1SShri Abhyankar Level: intermediate 339*5f2c45f1SShri Abhyankar 340*5f2c45f1SShri Abhyankar .seealso: DMNetworkGetNumComponents, DMNetworkGetComponentDataArray, 341*5f2c45f1SShri Abhyankar @*/ 342*5f2c45f1SShri Abhyankar PetscErrorCode DMNetworkGetComponentTypeOffset(DM dm,PetscInt p, PetscInt compnum, PetscInt *compkey, PetscInt *offset) 343*5f2c45f1SShri Abhyankar { 344*5f2c45f1SShri Abhyankar PetscErrorCode ierr; 345*5f2c45f1SShri Abhyankar PetscInt offsetp; 346*5f2c45f1SShri Abhyankar DMNetworkComponentHeader header; 347*5f2c45f1SShri Abhyankar DM_Network *network = (DM_Network*)dm->data; 348*5f2c45f1SShri Abhyankar 349*5f2c45f1SShri Abhyankar PetscFunctionBegin; 350*5f2c45f1SShri Abhyankar ierr = PetscSectionGetOffset(network->DataSection,p,&offsetp);CHKERRQ(ierr); 351*5f2c45f1SShri Abhyankar header = (DMNetworkComponentHeader)(network->componentdataarray+offsetp); 352*5f2c45f1SShri Abhyankar *compkey = header->key[compnum]; 353*5f2c45f1SShri Abhyankar *offset = offsetp+network->dataheadersize+header->offset[compnum]; 354*5f2c45f1SShri Abhyankar PetscFunctionReturn(0); 355*5f2c45f1SShri Abhyankar } 356*5f2c45f1SShri Abhyankar 357*5f2c45f1SShri Abhyankar #undef __FUNCT__ 358*5f2c45f1SShri Abhyankar #define __FUNCT__ "DMNetworkGetVariableOffset" 359*5f2c45f1SShri Abhyankar /*@ 360*5f2c45f1SShri Abhyankar DMNetworkGetVariableOffset - Get the offset for accessing the variable associated with the given vertex/edge from the local vector. 361*5f2c45f1SShri Abhyankar 362*5f2c45f1SShri Abhyankar Not Collective 363*5f2c45f1SShri Abhyankar 364*5f2c45f1SShri Abhyankar Input Parameters: 365*5f2c45f1SShri Abhyankar + dm - The DMNetwork object 366*5f2c45f1SShri Abhyankar - p - the edge/vertex point 367*5f2c45f1SShri Abhyankar 368*5f2c45f1SShri Abhyankar Output Parameters: 369*5f2c45f1SShri Abhyankar . offset - the offset 370*5f2c45f1SShri Abhyankar 371*5f2c45f1SShri Abhyankar Level: intermediate 372*5f2c45f1SShri Abhyankar 373*5f2c45f1SShri Abhyankar .seealso: DMNetworkGetVariableGlobalOffset, DMGetLocalVector 374*5f2c45f1SShri Abhyankar @*/ 375*5f2c45f1SShri Abhyankar PetscErrorCode DMNetworkGetVariableOffset(DM dm,PetscInt p,PetscInt *offset) 376*5f2c45f1SShri Abhyankar { 377*5f2c45f1SShri Abhyankar PetscErrorCode ierr; 378*5f2c45f1SShri Abhyankar DM_Network *network = (DM_Network*)dm->data; 379*5f2c45f1SShri Abhyankar 380*5f2c45f1SShri Abhyankar PetscFunctionBegin; 381*5f2c45f1SShri Abhyankar ierr = PetscSectionGetOffset(network->DofSection,p,offset);CHKERRQ(ierr); 382*5f2c45f1SShri Abhyankar PetscFunctionReturn(0); 383*5f2c45f1SShri Abhyankar } 384*5f2c45f1SShri Abhyankar 385*5f2c45f1SShri Abhyankar #undef __FUNCT__ 386*5f2c45f1SShri Abhyankar #define __FUNCT__ "DMNetworkGetVariableGlobalOffset" 387*5f2c45f1SShri Abhyankar /*@ 388*5f2c45f1SShri Abhyankar DMNetworkGetVariableGlobalOffset - Get the global offset for the variable associated with the given vertex/edge from the global vector. 389*5f2c45f1SShri Abhyankar 390*5f2c45f1SShri Abhyankar Not Collective 391*5f2c45f1SShri Abhyankar 392*5f2c45f1SShri Abhyankar Input Parameters: 393*5f2c45f1SShri Abhyankar + dm - The DMNetwork object 394*5f2c45f1SShri Abhyankar - p - the edge/vertex point 395*5f2c45f1SShri Abhyankar 396*5f2c45f1SShri Abhyankar Output Parameters: 397*5f2c45f1SShri Abhyankar . offsetg - the offset 398*5f2c45f1SShri Abhyankar 399*5f2c45f1SShri Abhyankar Level: intermediate 400*5f2c45f1SShri Abhyankar 401*5f2c45f1SShri Abhyankar .seealso: DMNetworkGetVariableOffset, DMGetLocalVector 402*5f2c45f1SShri Abhyankar @*/ 403*5f2c45f1SShri Abhyankar PetscErrorCode DMNetworkGetVariableGlobalOffset(DM dm,PetscInt p,PetscInt *offsetg) 404*5f2c45f1SShri Abhyankar { 405*5f2c45f1SShri Abhyankar PetscErrorCode ierr; 406*5f2c45f1SShri Abhyankar DM_Network *network = (DM_Network*)dm->data; 407*5f2c45f1SShri Abhyankar 408*5f2c45f1SShri Abhyankar PetscFunctionBegin; 409*5f2c45f1SShri Abhyankar ierr = PetscSectionGetOffset(network->GlobalDofSection,p,offsetg);CHKERRQ(ierr); 410*5f2c45f1SShri Abhyankar PetscFunctionReturn(0); 411*5f2c45f1SShri Abhyankar } 412*5f2c45f1SShri Abhyankar 413*5f2c45f1SShri Abhyankar #undef __FUNCT__ 414*5f2c45f1SShri Abhyankar #define __FUNCT__ "DMNetworkAddNumVariables" 415*5f2c45f1SShri Abhyankar /*@ 416*5f2c45f1SShri Abhyankar DMNetworkAddNumVariables - Add number of variables associated with a given point. 417*5f2c45f1SShri Abhyankar 418*5f2c45f1SShri Abhyankar Not Collective 419*5f2c45f1SShri Abhyankar 420*5f2c45f1SShri Abhyankar Input Parameters: 421*5f2c45f1SShri Abhyankar + dm - The DMNetworkObject 422*5f2c45f1SShri Abhyankar . p - the vertex/edge point 423*5f2c45f1SShri Abhyankar - nvar - number of additional variables 424*5f2c45f1SShri Abhyankar 425*5f2c45f1SShri Abhyankar Level: intermediate 426*5f2c45f1SShri Abhyankar 427*5f2c45f1SShri Abhyankar .seealso: DMNetworkSetNumVariables 428*5f2c45f1SShri Abhyankar @*/ 429*5f2c45f1SShri Abhyankar PetscErrorCode DMNetworkAddNumVariables(DM dm,PetscInt p,PetscInt nvar) 430*5f2c45f1SShri Abhyankar { 431*5f2c45f1SShri Abhyankar PetscErrorCode ierr; 432*5f2c45f1SShri Abhyankar DM_Network *network = (DM_Network*)dm->data; 433*5f2c45f1SShri Abhyankar 434*5f2c45f1SShri Abhyankar PetscFunctionBegin; 435*5f2c45f1SShri Abhyankar ierr = PetscSectionAddDof(network->DofSection,p,nvar);CHKERRQ(ierr); 436*5f2c45f1SShri Abhyankar PetscFunctionReturn(0); 437*5f2c45f1SShri Abhyankar } 438*5f2c45f1SShri Abhyankar 439*5f2c45f1SShri Abhyankar #undef __FUNCT__ 440*5f2c45f1SShri Abhyankar #define __FUNCT__ "DMNetworkAddNumVariables" 441*5f2c45f1SShri Abhyankar /*@ 442*5f2c45f1SShri Abhyankar DMNetworkSetNumVariables - Sets number of variables for a vertex/edge point. 443*5f2c45f1SShri Abhyankar 444*5f2c45f1SShri Abhyankar Not Collective 445*5f2c45f1SShri Abhyankar 446*5f2c45f1SShri Abhyankar Input Parameters: 447*5f2c45f1SShri Abhyankar + dm - The DMNetworkObject 448*5f2c45f1SShri Abhyankar . p - the vertex/edge point 449*5f2c45f1SShri Abhyankar - nvar - number of variables 450*5f2c45f1SShri Abhyankar 451*5f2c45f1SShri Abhyankar Level: intermediate 452*5f2c45f1SShri Abhyankar 453*5f2c45f1SShri Abhyankar .seealso: DMNetworkAddNumVariables 454*5f2c45f1SShri Abhyankar @*/ 455*5f2c45f1SShri Abhyankar PetscErrorCode DMNetworkSetNumVariables(DM dm,PetscInt p,PetscInt nvar) 456*5f2c45f1SShri Abhyankar { 457*5f2c45f1SShri Abhyankar PetscErrorCode ierr; 458*5f2c45f1SShri Abhyankar DM_Network *network = (DM_Network*)dm->data; 459*5f2c45f1SShri Abhyankar 460*5f2c45f1SShri Abhyankar PetscFunctionBegin; 461*5f2c45f1SShri Abhyankar ierr = PetscSectionSetDof(network->DofSection,p,nvar);CHKERRQ(ierr); 462*5f2c45f1SShri Abhyankar PetscFunctionReturn(0); 463*5f2c45f1SShri Abhyankar } 464*5f2c45f1SShri Abhyankar 465*5f2c45f1SShri Abhyankar /* Sets up the array that holds the data for all components and its associated section. This 466*5f2c45f1SShri Abhyankar function is called during DMSetUp() */ 467*5f2c45f1SShri Abhyankar #undef __FUNCT__ 468*5f2c45f1SShri Abhyankar #define __FUNCT__ "DMNetworkComponentSetUp" 469*5f2c45f1SShri Abhyankar PetscErrorCode DMNetworkComponentSetUp(DM dm) 470*5f2c45f1SShri Abhyankar { 471*5f2c45f1SShri Abhyankar PetscErrorCode ierr; 472*5f2c45f1SShri Abhyankar DM_Network *network = (DM_Network*)dm->data; 473*5f2c45f1SShri Abhyankar PetscInt arr_size; 474*5f2c45f1SShri Abhyankar PetscInt p,offset,offsetp; 475*5f2c45f1SShri Abhyankar DMNetworkComponentHeader header; 476*5f2c45f1SShri Abhyankar DMNetworkComponentValue cvalue; 477*5f2c45f1SShri Abhyankar DMNetworkComponentGenericDataType *componentdataarray; 478*5f2c45f1SShri Abhyankar PetscInt ncomp, i; 479*5f2c45f1SShri Abhyankar 480*5f2c45f1SShri Abhyankar PetscFunctionBegin; 481*5f2c45f1SShri Abhyankar ierr = PetscSectionSetUp(network->DataSection);CHKERRQ(ierr); 482*5f2c45f1SShri Abhyankar ierr = PetscSectionGetStorageSize(network->DataSection,&arr_size);CHKERRQ(ierr); 483*5f2c45f1SShri Abhyankar ierr = PetscMalloc(arr_size*sizeof(DMNetworkComponentGenericDataType),&network->componentdataarray);CHKERRQ(ierr); 484*5f2c45f1SShri Abhyankar componentdataarray = network->componentdataarray; 485*5f2c45f1SShri Abhyankar for (p = network->pStart; p < network->pEnd; p++) { 486*5f2c45f1SShri Abhyankar ierr = PetscSectionGetOffset(network->DataSection,p,&offsetp);CHKERRQ(ierr); 487*5f2c45f1SShri Abhyankar /* Copy header */ 488*5f2c45f1SShri Abhyankar header = &network->header[p]; 489*5f2c45f1SShri Abhyankar ierr = PetscMemcpy(componentdataarray+offsetp,header,network->dataheadersize*sizeof(DMNetworkComponentGenericDataType)); 490*5f2c45f1SShri Abhyankar /* Copy data */ 491*5f2c45f1SShri Abhyankar cvalue = &network->cvalue[p]; 492*5f2c45f1SShri Abhyankar ncomp = header->ndata; 493*5f2c45f1SShri Abhyankar for (i = 0; i < ncomp; i++) { 494*5f2c45f1SShri Abhyankar offset = offsetp + network->dataheadersize + header->offset[i]; 495*5f2c45f1SShri Abhyankar ierr = PetscMemcpy(componentdataarray+offset,cvalue->data[i],header->size[i]*sizeof(DMNetworkComponentGenericDataType)); 496*5f2c45f1SShri Abhyankar } 497*5f2c45f1SShri Abhyankar } 498*5f2c45f1SShri Abhyankar PetscFunctionReturn(0); 499*5f2c45f1SShri Abhyankar } 500*5f2c45f1SShri Abhyankar 501*5f2c45f1SShri Abhyankar /* Sets up the section for dofs. This routine is called during DMSetUp() */ 502*5f2c45f1SShri Abhyankar #undef __FUNCT__ 503*5f2c45f1SShri Abhyankar #define __FUNCT__ "DMNetworkVariablesSetUp" 504*5f2c45f1SShri Abhyankar PetscErrorCode DMNetworkVariablesSetUp(DM dm) 505*5f2c45f1SShri Abhyankar { 506*5f2c45f1SShri Abhyankar PetscErrorCode ierr; 507*5f2c45f1SShri Abhyankar DM_Network *network = (DM_Network*)dm->data; 508*5f2c45f1SShri Abhyankar 509*5f2c45f1SShri Abhyankar PetscFunctionBegin; 510*5f2c45f1SShri Abhyankar ierr = PetscSectionSetUp(network->DofSection);CHKERRQ(ierr); 511*5f2c45f1SShri Abhyankar PetscFunctionReturn(0); 512*5f2c45f1SShri Abhyankar } 513*5f2c45f1SShri Abhyankar 514*5f2c45f1SShri Abhyankar #undef __FUNCT__ 515*5f2c45f1SShri Abhyankar #define __FUNCT__ "DMNetworkGetComponentDataArray" 516*5f2c45f1SShri Abhyankar /*@C 517*5f2c45f1SShri Abhyankar DMNetworkGetComponentDataArray - Returns the component data array 518*5f2c45f1SShri Abhyankar 519*5f2c45f1SShri Abhyankar Not Collective 520*5f2c45f1SShri Abhyankar 521*5f2c45f1SShri Abhyankar Input Parameters: 522*5f2c45f1SShri Abhyankar . dm - The DMNetwork Object 523*5f2c45f1SShri Abhyankar 524*5f2c45f1SShri Abhyankar Output Parameters: 525*5f2c45f1SShri Abhyankar . componentdataarray - array that holds data for all components 526*5f2c45f1SShri Abhyankar 527*5f2c45f1SShri Abhyankar Level: intermediate 528*5f2c45f1SShri Abhyankar 529*5f2c45f1SShri Abhyankar .seealso: DMNetworkGetComponentTypeOffset, DMNetworkGetNumComponents 530*5f2c45f1SShri Abhyankar @*/ 531*5f2c45f1SShri Abhyankar PetscErrorCode DMNetworkGetComponentDataArray(DM dm,DMNetworkComponentGenericDataType **componentdataarray) 532*5f2c45f1SShri Abhyankar { 533*5f2c45f1SShri Abhyankar DM_Network *network = (DM_Network*)dm->data; 534*5f2c45f1SShri Abhyankar 535*5f2c45f1SShri Abhyankar PetscFunctionBegin; 536*5f2c45f1SShri Abhyankar *componentdataarray = network->componentdataarray; 537*5f2c45f1SShri Abhyankar PetscFunctionReturn(0); 538*5f2c45f1SShri Abhyankar } 539*5f2c45f1SShri Abhyankar 540*5f2c45f1SShri Abhyankar #undef __FUNCT__ 541*5f2c45f1SShri Abhyankar #define __FUNCT__ "DMNetworkDistribute" 542*5f2c45f1SShri Abhyankar /*@ 543*5f2c45f1SShri Abhyankar DMNetworkDistribute - Distributes the network and moves associated component data. 544*5f2c45f1SShri Abhyankar 545*5f2c45f1SShri Abhyankar Collective 546*5f2c45f1SShri Abhyankar 547*5f2c45f1SShri Abhyankar Input Parameter: 548*5f2c45f1SShri Abhyankar + oldDM - the original DMNetwork object 549*5f2c45f1SShri Abhyankar . partitioner - The partitioning package, or NULL for the default 550*5f2c45f1SShri Abhyankar - overlap - The overlap of partitions, 0 is the default 551*5f2c45f1SShri Abhyankar 552*5f2c45f1SShri Abhyankar Output Parameter: 553*5f2c45f1SShri Abhyankar . distDM - the distributed DMNetwork object 554*5f2c45f1SShri Abhyankar 555*5f2c45f1SShri Abhyankar Notes: 556*5f2c45f1SShri Abhyankar This routine should be called only when using multiple processors. 557*5f2c45f1SShri Abhyankar 558*5f2c45f1SShri Abhyankar Distributes the network with a non-overlapping partitioning of the edges. 559*5f2c45f1SShri Abhyankar 560*5f2c45f1SShri Abhyankar Level: intermediate 561*5f2c45f1SShri Abhyankar 562*5f2c45f1SShri Abhyankar .seealso: DMNetworkCreate 563*5f2c45f1SShri Abhyankar @*/ 564*5f2c45f1SShri Abhyankar PetscErrorCode DMNetworkDistribute(DM oldDM, const char partitioner[], PetscInt overlap,DM *distDM) 565*5f2c45f1SShri Abhyankar { 566*5f2c45f1SShri Abhyankar PetscErrorCode ierr; 567*5f2c45f1SShri Abhyankar DM_Network *oldDMnetwork = (DM_Network*)oldDM->data; 568*5f2c45f1SShri Abhyankar PetscSF pointsf; 569*5f2c45f1SShri Abhyankar DM newDM; 570*5f2c45f1SShri Abhyankar DM_Network *newDMnetwork; 571*5f2c45f1SShri Abhyankar 572*5f2c45f1SShri Abhyankar PetscFunctionBegin; 573*5f2c45f1SShri Abhyankar ierr = DMNetworkCreate(PetscObjectComm((PetscObject)oldDM),&newDM);CHKERRQ(ierr); 574*5f2c45f1SShri Abhyankar newDMnetwork = (DM_Network*)newDM->data; 575*5f2c45f1SShri Abhyankar newDMnetwork->dataheadersize = sizeof(struct _p_DMNetworkComponentHeader)/sizeof(DMNetworkComponentGenericDataType); 576*5f2c45f1SShri Abhyankar /* Distribute plex dm and dof section */ 577*5f2c45f1SShri Abhyankar ierr = DMPlexDistribute(oldDMnetwork->plex,partitioner,overlap,&pointsf,&newDMnetwork->plex);CHKERRQ(ierr); 578*5f2c45f1SShri Abhyankar /* Distribute dof section */ 579*5f2c45f1SShri Abhyankar ierr = PetscSectionCreate(PetscObjectComm((PetscObject)oldDM),&newDMnetwork->DofSection);CHKERRQ(ierr); 580*5f2c45f1SShri Abhyankar ierr = PetscSFDistributeSection(pointsf,oldDMnetwork->DofSection,NULL,newDMnetwork->DofSection);CHKERRQ(ierr); 581*5f2c45f1SShri Abhyankar ierr = PetscSectionCreate(PetscObjectComm((PetscObject)oldDM),&newDMnetwork->DataSection);CHKERRQ(ierr); 582*5f2c45f1SShri Abhyankar /* Distribute data and associated section */ 583*5f2c45f1SShri Abhyankar ierr = DMPlexDistributeData(newDMnetwork->plex,pointsf,oldDMnetwork->DataSection,MPI_INT,(void*)oldDMnetwork->componentdataarray,newDMnetwork->DataSection,(void**)&newDMnetwork->componentdataarray);CHKERRQ(ierr); 584*5f2c45f1SShri Abhyankar /* Destroy point SF */ 585*5f2c45f1SShri Abhyankar ierr = PetscSFDestroy(&pointsf);CHKERRQ(ierr); 586*5f2c45f1SShri Abhyankar 587*5f2c45f1SShri Abhyankar ierr = PetscSectionGetChart(newDMnetwork->DataSection,&newDMnetwork->pStart,&newDMnetwork->pEnd);CHKERRQ(ierr); 588*5f2c45f1SShri Abhyankar ierr = DMPlexGetHeightStratum(newDMnetwork->plex,0, &newDMnetwork->eStart,&newDMnetwork->eEnd);CHKERRQ(ierr); 589*5f2c45f1SShri Abhyankar ierr = DMPlexGetHeightStratum(newDMnetwork->plex,1,&newDMnetwork->vStart,&newDMnetwork->vEnd);CHKERRQ(ierr); 590*5f2c45f1SShri Abhyankar newDMnetwork->nEdges = newDMnetwork->eEnd - newDMnetwork->eStart; 591*5f2c45f1SShri Abhyankar newDMnetwork->nNodes = newDMnetwork->vEnd - newDMnetwork->vStart; 592*5f2c45f1SShri Abhyankar newDMnetwork->NNodes = oldDMnetwork->NNodes; 593*5f2c45f1SShri Abhyankar newDMnetwork->NEdges = oldDMnetwork->NEdges; 594*5f2c45f1SShri Abhyankar /* Set Dof section as the default section for dm */ 595*5f2c45f1SShri Abhyankar ierr = DMSetDefaultSection(newDMnetwork->plex,newDMnetwork->DofSection);CHKERRQ(ierr); 596*5f2c45f1SShri Abhyankar ierr = DMGetDefaultGlobalSection(newDMnetwork->plex,&newDMnetwork->GlobalDofSection);CHKERRQ(ierr); 597*5f2c45f1SShri Abhyankar 598*5f2c45f1SShri Abhyankar *distDM = newDM; 599*5f2c45f1SShri Abhyankar PetscFunctionReturn(0); 600*5f2c45f1SShri Abhyankar } 601*5f2c45f1SShri Abhyankar 602*5f2c45f1SShri Abhyankar #undef __FUNCT__ 603*5f2c45f1SShri Abhyankar #define __FUNCT__ "DMNetworkGetSupportingEdges" 604*5f2c45f1SShri Abhyankar /*@C 605*5f2c45f1SShri Abhyankar DMNetworkGetSupportingEdges - Return the supporting edges for this vertex point 606*5f2c45f1SShri Abhyankar 607*5f2c45f1SShri Abhyankar Not Collective 608*5f2c45f1SShri Abhyankar 609*5f2c45f1SShri Abhyankar Input Parameters: 610*5f2c45f1SShri Abhyankar + dm - The DMNetwork object 611*5f2c45f1SShri Abhyankar - p - the vertex point 612*5f2c45f1SShri Abhyankar 613*5f2c45f1SShri Abhyankar Output Paramters: 614*5f2c45f1SShri Abhyankar + nedges - number of edges connected to this vertex point 615*5f2c45f1SShri Abhyankar - edges - List of edge points 616*5f2c45f1SShri Abhyankar 617*5f2c45f1SShri Abhyankar Level: intermediate 618*5f2c45f1SShri Abhyankar 619*5f2c45f1SShri Abhyankar Fortran Notes: 620*5f2c45f1SShri Abhyankar Since it returns an array, this routine is only available in Fortran 90, and you must 621*5f2c45f1SShri Abhyankar include petsc.h90 in your code. 622*5f2c45f1SShri Abhyankar 623*5f2c45f1SShri Abhyankar .seealso: DMNetworkCreate, DMNetworkGetConnectedNodes 624*5f2c45f1SShri Abhyankar @*/ 625*5f2c45f1SShri Abhyankar PetscErrorCode DMNetworkGetSupportingEdges(DM dm,PetscInt vertex,PetscInt *nedges,const PetscInt *edges[]) 626*5f2c45f1SShri Abhyankar { 627*5f2c45f1SShri Abhyankar PetscErrorCode ierr; 628*5f2c45f1SShri Abhyankar DM_Network *network = (DM_Network*)dm->data; 629*5f2c45f1SShri Abhyankar 630*5f2c45f1SShri Abhyankar PetscFunctionBegin; 631*5f2c45f1SShri Abhyankar ierr = DMPlexGetSupportSize(network->plex,vertex,nedges);CHKERRQ(ierr); 632*5f2c45f1SShri Abhyankar ierr = DMPlexGetSupport(network->plex,vertex,edges);CHKERRQ(ierr); 633*5f2c45f1SShri Abhyankar PetscFunctionReturn(0); 634*5f2c45f1SShri Abhyankar } 635*5f2c45f1SShri Abhyankar 636*5f2c45f1SShri Abhyankar #undef __FUNCT__ 637*5f2c45f1SShri Abhyankar #define __FUNCT__ "DMNetworkGetConnectedNodes" 638*5f2c45f1SShri Abhyankar /*@C 639*5f2c45f1SShri Abhyankar DMNetworkGetConnectedNodes - Return the connected edges for this edge point 640*5f2c45f1SShri Abhyankar 641*5f2c45f1SShri Abhyankar Not Collective 642*5f2c45f1SShri Abhyankar 643*5f2c45f1SShri Abhyankar Input Parameters: 644*5f2c45f1SShri Abhyankar + dm - The DMNetwork object 645*5f2c45f1SShri Abhyankar - p - the edge point 646*5f2c45f1SShri Abhyankar 647*5f2c45f1SShri Abhyankar Output Paramters: 648*5f2c45f1SShri Abhyankar . vertices - vertices connected to this edge 649*5f2c45f1SShri Abhyankar 650*5f2c45f1SShri Abhyankar Level: intermediate 651*5f2c45f1SShri Abhyankar 652*5f2c45f1SShri Abhyankar Fortran Notes: 653*5f2c45f1SShri Abhyankar Since it returns an array, this routine is only available in Fortran 90, and you must 654*5f2c45f1SShri Abhyankar include petsc.h90 in your code. 655*5f2c45f1SShri Abhyankar 656*5f2c45f1SShri Abhyankar .seealso: DMNetworkCreate, DMNetworkGetSupportingEdges 657*5f2c45f1SShri Abhyankar @*/ 658*5f2c45f1SShri Abhyankar PetscErrorCode DMNetworkGetConnectedNodes(DM dm,PetscInt edge,const PetscInt *vertices[]) 659*5f2c45f1SShri Abhyankar { 660*5f2c45f1SShri Abhyankar PetscErrorCode ierr; 661*5f2c45f1SShri Abhyankar DM_Network *network = (DM_Network*)dm->data; 662*5f2c45f1SShri Abhyankar 663*5f2c45f1SShri Abhyankar PetscFunctionBegin; 664*5f2c45f1SShri Abhyankar ierr = DMPlexGetCone(network->plex,edge,vertices);CHKERRQ(ierr); 665*5f2c45f1SShri Abhyankar PetscFunctionReturn(0); 666*5f2c45f1SShri Abhyankar } 667*5f2c45f1SShri Abhyankar 668*5f2c45f1SShri Abhyankar #undef __FUNCT__ 669*5f2c45f1SShri Abhyankar #define __FUNCT__ "DMNetworkIsGhostVertex" 670*5f2c45f1SShri Abhyankar /*@ 671*5f2c45f1SShri Abhyankar DMNetworkIsGhostVertex - Returns TRUE if the vertex is a ghost vertex 672*5f2c45f1SShri Abhyankar 673*5f2c45f1SShri Abhyankar Not Collective 674*5f2c45f1SShri Abhyankar 675*5f2c45f1SShri Abhyankar Input Parameters: 676*5f2c45f1SShri Abhyankar + dm - The DMNetwork object 677*5f2c45f1SShri Abhyankar . p - the vertex point 678*5f2c45f1SShri Abhyankar 679*5f2c45f1SShri Abhyankar Output Parameter: 680*5f2c45f1SShri Abhyankar . isghost - TRUE if the vertex is a ghost point 681*5f2c45f1SShri Abhyankar 682*5f2c45f1SShri Abhyankar Level: intermediate 683*5f2c45f1SShri Abhyankar 684*5f2c45f1SShri Abhyankar .seealso: DMNetworkCreate, DMNetworkGetConnectedNodes, DMNetworkGetVertexRange 685*5f2c45f1SShri Abhyankar @*/ 686*5f2c45f1SShri Abhyankar PetscErrorCode DMNetworkIsGhostVertex(DM dm,PetscInt p,PetscBool *isghost) 687*5f2c45f1SShri Abhyankar { 688*5f2c45f1SShri Abhyankar PetscErrorCode ierr; 689*5f2c45f1SShri Abhyankar DM_Network *network = (DM_Network*)dm->data; 690*5f2c45f1SShri Abhyankar PetscInt offsetg; 691*5f2c45f1SShri Abhyankar PetscSection sectiong; 692*5f2c45f1SShri Abhyankar 693*5f2c45f1SShri Abhyankar PetscFunctionBegin; 694*5f2c45f1SShri Abhyankar *isghost = PETSC_FALSE; 695*5f2c45f1SShri Abhyankar ierr = DMGetDefaultGlobalSection(network->plex,§iong);CHKERRQ(ierr); 696*5f2c45f1SShri Abhyankar ierr = PetscSectionGetOffset(sectiong,p,&offsetg);CHKERRQ(ierr); 697*5f2c45f1SShri Abhyankar if (offsetg < 0) *isghost = PETSC_TRUE; 698*5f2c45f1SShri Abhyankar PetscFunctionReturn(0); 699*5f2c45f1SShri Abhyankar } 700*5f2c45f1SShri Abhyankar 701*5f2c45f1SShri Abhyankar #undef __FUNCT__ 702*5f2c45f1SShri Abhyankar #define __FUNCT__ "DMSetUp_Network" 703*5f2c45f1SShri Abhyankar PetscErrorCode DMSetUp_Network(DM dm) 704*5f2c45f1SShri Abhyankar { 705*5f2c45f1SShri Abhyankar PetscErrorCode ierr; 706*5f2c45f1SShri Abhyankar DM_Network *network=(DM_Network*)dm->data; 707*5f2c45f1SShri Abhyankar 708*5f2c45f1SShri Abhyankar PetscFunctionBegin; 709*5f2c45f1SShri Abhyankar ierr = DMNetworkComponentSetUp(dm);CHKERRQ(ierr); 710*5f2c45f1SShri Abhyankar ierr = DMNetworkVariablesSetUp(dm);CHKERRQ(ierr); 711*5f2c45f1SShri Abhyankar 712*5f2c45f1SShri Abhyankar ierr = DMSetDefaultSection(network->plex,network->DofSection);CHKERRQ(ierr); 713*5f2c45f1SShri Abhyankar ierr = DMGetDefaultGlobalSection(network->plex,&network->GlobalDofSection);CHKERRQ(ierr); 714*5f2c45f1SShri Abhyankar PetscFunctionReturn(0); 715*5f2c45f1SShri Abhyankar } 716*5f2c45f1SShri Abhyankar 717*5f2c45f1SShri Abhyankar #undef __FUNCT__ 718*5f2c45f1SShri Abhyankar #define __FUNCT__ "DMCreateMatrix_Network" 719*5f2c45f1SShri Abhyankar PetscErrorCode DMCreateMatrix_Network(DM dm,Mat *J) 720*5f2c45f1SShri Abhyankar { 721*5f2c45f1SShri Abhyankar PetscErrorCode ierr; 722*5f2c45f1SShri Abhyankar DM_Network *network = (DM_Network*) dm->data; 723*5f2c45f1SShri Abhyankar 724*5f2c45f1SShri Abhyankar PetscFunctionBegin; 725*5f2c45f1SShri Abhyankar ierr = DMCreateMatrix(network->plex,J);CHKERRQ(ierr); 726*5f2c45f1SShri Abhyankar ierr = MatSetDM(*J,dm);CHKERRQ(ierr); 727*5f2c45f1SShri Abhyankar PetscFunctionReturn(0); 728*5f2c45f1SShri Abhyankar } 729*5f2c45f1SShri Abhyankar 730*5f2c45f1SShri Abhyankar #undef __FUNCT__ 731*5f2c45f1SShri Abhyankar #define __FUNCT__ "DMDestroy_Network" 732*5f2c45f1SShri Abhyankar PetscErrorCode DMDestroy_Network(DM dm) 733*5f2c45f1SShri Abhyankar { 734*5f2c45f1SShri Abhyankar PetscErrorCode ierr; 735*5f2c45f1SShri Abhyankar DM_Network *network = (DM_Network*) dm->data; 736*5f2c45f1SShri Abhyankar 737*5f2c45f1SShri Abhyankar PetscFunctionBegin; 738*5f2c45f1SShri Abhyankar ierr = DMDestroy(&network->plex);CHKERRQ(ierr); 739*5f2c45f1SShri Abhyankar network->edges = NULL; 740*5f2c45f1SShri Abhyankar ierr = PetscSectionDestroy(&network->DataSection);CHKERRQ(ierr); 741*5f2c45f1SShri Abhyankar ierr = PetscSectionDestroy(&network->DofSection);CHKERRQ(ierr); 742*5f2c45f1SShri Abhyankar /* ierr = PetscSectionDestroy(&network->GlobalDofSection);CHKERRQ(ierr); */ 743*5f2c45f1SShri Abhyankar ierr = PetscFree(network->componentdataarray);CHKERRQ(ierr); 744*5f2c45f1SShri Abhyankar ierr = PetscFree(network->cvalue);CHKERRQ(ierr); 745*5f2c45f1SShri Abhyankar ierr = PetscFree(network->header);CHKERRQ(ierr); 746*5f2c45f1SShri Abhyankar ierr = PetscFree(network);CHKERRQ(ierr); 747*5f2c45f1SShri Abhyankar PetscFunctionReturn(0); 748*5f2c45f1SShri Abhyankar } 749*5f2c45f1SShri Abhyankar 750*5f2c45f1SShri Abhyankar #undef __FUNCT__ 751*5f2c45f1SShri Abhyankar #define __FUNCT__ "DMView_Network" 752*5f2c45f1SShri Abhyankar PetscErrorCode DMView_Network(DM dm, PetscViewer viewer) 753*5f2c45f1SShri Abhyankar { 754*5f2c45f1SShri Abhyankar PetscErrorCode ierr; 755*5f2c45f1SShri Abhyankar DM_Network *network = (DM_Network*) dm->data; 756*5f2c45f1SShri Abhyankar 757*5f2c45f1SShri Abhyankar PetscFunctionBegin; 758*5f2c45f1SShri Abhyankar ierr = DMView(network->plex,viewer);CHKERRQ(ierr); 759*5f2c45f1SShri Abhyankar PetscFunctionReturn(0); 760*5f2c45f1SShri Abhyankar } 761*5f2c45f1SShri Abhyankar 762*5f2c45f1SShri Abhyankar #undef __FUNCT__ 763*5f2c45f1SShri Abhyankar #define __FUNCT__ "DMGlobalToLocalBegin_Network" 764*5f2c45f1SShri Abhyankar PetscErrorCode DMGlobalToLocalBegin_Network(DM dm, Vec g, InsertMode mode, Vec l) 765*5f2c45f1SShri Abhyankar { 766*5f2c45f1SShri Abhyankar PetscErrorCode ierr; 767*5f2c45f1SShri Abhyankar DM_Network *network = (DM_Network*) dm->data; 768*5f2c45f1SShri Abhyankar 769*5f2c45f1SShri Abhyankar PetscFunctionBegin; 770*5f2c45f1SShri Abhyankar ierr = DMGlobalToLocalBegin(network->plex,g,mode,l);CHKERRQ(ierr); 771*5f2c45f1SShri Abhyankar PetscFunctionReturn(0); 772*5f2c45f1SShri Abhyankar } 773*5f2c45f1SShri Abhyankar 774*5f2c45f1SShri Abhyankar #undef __FUNCT__ 775*5f2c45f1SShri Abhyankar #define __FUNCT__ "DMGlobalToLocalEnd_Network" 776*5f2c45f1SShri Abhyankar PetscErrorCode DMGlobalToLocalEnd_Network(DM dm, Vec g, InsertMode mode, Vec l) 777*5f2c45f1SShri Abhyankar { 778*5f2c45f1SShri Abhyankar PetscErrorCode ierr; 779*5f2c45f1SShri Abhyankar DM_Network *network = (DM_Network*) dm->data; 780*5f2c45f1SShri Abhyankar 781*5f2c45f1SShri Abhyankar PetscFunctionBegin; 782*5f2c45f1SShri Abhyankar ierr = DMGlobalToLocalEnd(network->plex,g,mode,l);CHKERRQ(ierr); 783*5f2c45f1SShri Abhyankar PetscFunctionReturn(0); 784*5f2c45f1SShri Abhyankar } 785*5f2c45f1SShri Abhyankar 786*5f2c45f1SShri Abhyankar #undef __FUNCT__ 787*5f2c45f1SShri Abhyankar #define __FUNCT__ "DMLocalToGlobalBegin_Network" 788*5f2c45f1SShri Abhyankar PetscErrorCode DMLocalToGlobalBegin_Network(DM dm, Vec l, InsertMode mode, Vec g) 789*5f2c45f1SShri Abhyankar { 790*5f2c45f1SShri Abhyankar PetscErrorCode ierr; 791*5f2c45f1SShri Abhyankar DM_Network *network = (DM_Network*) dm->data; 792*5f2c45f1SShri Abhyankar 793*5f2c45f1SShri Abhyankar PetscFunctionBegin; 794*5f2c45f1SShri Abhyankar ierr = DMLocalToGlobalBegin(network->plex,l,mode,g);CHKERRQ(ierr); 795*5f2c45f1SShri Abhyankar PetscFunctionReturn(0); 796*5f2c45f1SShri Abhyankar } 797*5f2c45f1SShri Abhyankar 798*5f2c45f1SShri Abhyankar #undef __FUNCT__ 799*5f2c45f1SShri Abhyankar #define __FUNCT__ "DMLocalToGlobalEnd_Network" 800*5f2c45f1SShri Abhyankar PetscErrorCode DMLocalToGlobalEnd_Network(DM dm, Vec l, InsertMode mode, Vec g) 801*5f2c45f1SShri Abhyankar { 802*5f2c45f1SShri Abhyankar PetscErrorCode ierr; 803*5f2c45f1SShri Abhyankar DM_Network *network = (DM_Network*) dm->data; 804*5f2c45f1SShri Abhyankar 805*5f2c45f1SShri Abhyankar PetscFunctionBegin; 806*5f2c45f1SShri Abhyankar ierr = DMLocalToGlobalEnd(network->plex,l,mode,g);CHKERRQ(ierr); 807*5f2c45f1SShri Abhyankar PetscFunctionReturn(0); 808*5f2c45f1SShri Abhyankar } 809