xref: /petsc/src/dm/impls/redundant/dmredundant.c (revision 73c9e547bc1063479ca31fa7bde242fcd7d19d63)
1af0996ceSBarry Smith #include <petsc/private/dmimpl.h>
28ac4e037SJed Brown #include <petscdmredundant.h>   /*I      "petscdmredundant.h" I*/
38ac4e037SJed Brown 
48ac4e037SJed Brown typedef struct  {
5907376e6SBarry Smith   PetscMPIInt rank;                /* owner */
68ac4e037SJed Brown   PetscInt    N;                   /* total number of dofs */
78ac4e037SJed Brown   PetscInt    n;                   /* owned number of dofs, n=N on owner, n=0 on non-owners */
88ac4e037SJed Brown } DM_Redundant;
98ac4e037SJed Brown 
10b412c318SBarry Smith static PetscErrorCode DMCreateMatrix_Redundant(DM dm,Mat *J)
118ac4e037SJed Brown {
128ac4e037SJed Brown   DM_Redundant           *red = (DM_Redundant*)dm->data;
138ac4e037SJed Brown   PetscErrorCode         ierr;
1445b6f7e9SBarry Smith   ISLocalToGlobalMapping ltog;
15d2e2b695SJed Brown   PetscInt               i,rstart,rend,*cols;
16d2e2b695SJed Brown   PetscScalar            *vals;
178ac4e037SJed Brown 
188ac4e037SJed Brown   PetscFunctionBegin;
19ce94432eSBarry Smith   ierr = MatCreate(PetscObjectComm((PetscObject)dm),J);CHKERRQ(ierr);
208ac4e037SJed Brown   ierr = MatSetSizes(*J,red->n,red->n,red->N,red->N);CHKERRQ(ierr);
21b412c318SBarry Smith   ierr = MatSetType(*J,dm->mattype);CHKERRQ(ierr);
220298fd71SBarry Smith   ierr = MatSeqAIJSetPreallocation(*J,red->n,NULL);CHKERRQ(ierr);
230298fd71SBarry Smith   ierr = MatSeqBAIJSetPreallocation(*J,1,red->n,NULL);CHKERRQ(ierr);
240298fd71SBarry Smith   ierr = MatMPIAIJSetPreallocation(*J,red->n,NULL,red->N-red->n,NULL);CHKERRQ(ierr);
250298fd71SBarry Smith   ierr = MatMPIBAIJSetPreallocation(*J,1,red->n,NULL,red->N-red->n,NULL);CHKERRQ(ierr);
268ac4e037SJed Brown 
278ac4e037SJed Brown   ierr = DMGetLocalToGlobalMapping(dm,&ltog);CHKERRQ(ierr);
288ac4e037SJed Brown   ierr = MatSetLocalToGlobalMapping(*J,ltog,ltog);CHKERRQ(ierr);
29d2e2b695SJed Brown 
30dcca6d9dSJed Brown   ierr = PetscMalloc2(red->N,&cols,red->N,&vals);CHKERRQ(ierr);
31d2e2b695SJed Brown   for (i=0; i<red->N; i++) {
32d2e2b695SJed Brown     cols[i] = i;
33d2e2b695SJed Brown     vals[i] = 0.0;
34d2e2b695SJed Brown   }
35d2e2b695SJed Brown   ierr = MatGetOwnershipRange(*J,&rstart,&rend);CHKERRQ(ierr);
36d2e2b695SJed Brown   for (i=rstart; i<rend; i++) {
37d2e2b695SJed Brown     ierr = MatSetValues(*J,1,&i,red->N,cols,vals,INSERT_VALUES);CHKERRQ(ierr);
38d2e2b695SJed Brown   }
39d2e2b695SJed Brown   ierr = PetscFree2(cols,vals);CHKERRQ(ierr);
40d2e2b695SJed Brown   ierr = MatAssemblyBegin(*J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
41d2e2b695SJed Brown   ierr = MatAssemblyEnd(*J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
428ac4e037SJed Brown   PetscFunctionReturn(0);
438ac4e037SJed Brown }
448ac4e037SJed Brown 
458ac4e037SJed Brown static PetscErrorCode DMDestroy_Redundant(DM dm)
468ac4e037SJed Brown {
478ac4e037SJed Brown   PetscErrorCode ierr;
488ac4e037SJed Brown 
498ac4e037SJed Brown   PetscFunctionBegin;
50bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)dm,"DMRedundantSetSize_C",NULL);CHKERRQ(ierr);
51bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)dm,"DMRedundantGetSize_C",NULL);CHKERRQ(ierr);
5236623e74SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)dm,"DMSetUpGLVisViewer_C",NULL);CHKERRQ(ierr);
53435a35e8SMatthew G Knepley   /* This was originally freed in DMDestroy(), but that prevents reference counting of backend objects */
54435a35e8SMatthew G Knepley   ierr = PetscFree(dm->data);CHKERRQ(ierr);
558ac4e037SJed Brown   PetscFunctionReturn(0);
568ac4e037SJed Brown }
578ac4e037SJed Brown 
588ac4e037SJed Brown static PetscErrorCode DMCreateGlobalVector_Redundant(DM dm,Vec *gvec)
598ac4e037SJed Brown {
608ac4e037SJed Brown   PetscErrorCode         ierr;
618ac4e037SJed Brown   DM_Redundant           *red = (DM_Redundant*)dm->data;
6245b6f7e9SBarry Smith   ISLocalToGlobalMapping ltog;
638ac4e037SJed Brown 
648ac4e037SJed Brown   PetscFunctionBegin;
658ac4e037SJed Brown   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
668ac4e037SJed Brown   PetscValidPointer(gvec,2);
678ac4e037SJed Brown   *gvec = 0;
68ce94432eSBarry Smith   ierr  = VecCreate(PetscObjectComm((PetscObject)dm),gvec);CHKERRQ(ierr);
698ac4e037SJed Brown   ierr  = VecSetSizes(*gvec,red->n,red->N);CHKERRQ(ierr);
708ac4e037SJed Brown   ierr  = VecSetType(*gvec,dm->vectype);CHKERRQ(ierr);
718ac4e037SJed Brown   ierr  = DMGetLocalToGlobalMapping(dm,&ltog);CHKERRQ(ierr);
728ac4e037SJed Brown   ierr  = VecSetLocalToGlobalMapping(*gvec,ltog);CHKERRQ(ierr);
73c688c046SMatthew G Knepley   ierr  = VecSetDM(*gvec,dm);CHKERRQ(ierr);
748ac4e037SJed Brown   PetscFunctionReturn(0);
758ac4e037SJed Brown }
768ac4e037SJed Brown 
778ac4e037SJed Brown static PetscErrorCode DMCreateLocalVector_Redundant(DM dm,Vec *lvec)
788ac4e037SJed Brown {
798ac4e037SJed Brown   PetscErrorCode ierr;
808ac4e037SJed Brown   DM_Redundant   *red = (DM_Redundant*)dm->data;
818ac4e037SJed Brown 
828ac4e037SJed Brown   PetscFunctionBegin;
838ac4e037SJed Brown   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
848ac4e037SJed Brown   PetscValidPointer(lvec,2);
858ac4e037SJed Brown   *lvec = 0;
868ac4e037SJed Brown   ierr  = VecCreate(PETSC_COMM_SELF,lvec);CHKERRQ(ierr);
878ac4e037SJed Brown   ierr  = VecSetSizes(*lvec,red->N,red->N);CHKERRQ(ierr);
888ac4e037SJed Brown   ierr  = VecSetType(*lvec,dm->vectype);CHKERRQ(ierr);
89c688c046SMatthew G Knepley   ierr  = VecSetDM(*lvec,dm);CHKERRQ(ierr);
908ac4e037SJed Brown   PetscFunctionReturn(0);
918ac4e037SJed Brown }
928ac4e037SJed Brown 
938ac4e037SJed Brown static PetscErrorCode DMLocalToGlobalBegin_Redundant(DM dm,Vec l,InsertMode imode,Vec g)
948ac4e037SJed Brown {
958ac4e037SJed Brown   PetscErrorCode    ierr;
968ac4e037SJed Brown   DM_Redundant      *red = (DM_Redundant*)dm->data;
978ac4e037SJed Brown   const PetscScalar *lv;
988ac4e037SJed Brown   PetscScalar       *gv;
998ac4e037SJed Brown   PetscMPIInt       rank;
1008ac4e037SJed Brown 
1018ac4e037SJed Brown   PetscFunctionBegin;
102ce94432eSBarry Smith   ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)dm),&rank);CHKERRQ(ierr);
1038ac4e037SJed Brown   ierr = VecGetArrayRead(l,&lv);CHKERRQ(ierr);
1048ac4e037SJed Brown   ierr = VecGetArray(g,&gv);CHKERRQ(ierr);
1058ac4e037SJed Brown   switch (imode) {
1068ac4e037SJed Brown   case ADD_VALUES:
1078ac4e037SJed Brown   case MAX_VALUES:
1088ac4e037SJed Brown   {
1098ac4e037SJed Brown     void        *source;
1108ac4e037SJed Brown     PetscScalar *buffer;
1118ac4e037SJed Brown     PetscInt    i;
1128ac4e037SJed Brown     if (rank == red->rank) {
1138ac4e037SJed Brown #if defined(PETSC_HAVE_MPI_IN_PLACE)
1148ac4e037SJed Brown       buffer = gv;
1158ac4e037SJed Brown       source = MPI_IN_PLACE;
1168ac4e037SJed Brown #else
117785e854fSJed Brown       ierr   = PetscMalloc1(red->N,&buffer);CHKERRQ(ierr);
1188ac4e037SJed Brown       source = buffer;
1198ac4e037SJed Brown #endif
1208ac4e037SJed Brown       if (imode == ADD_VALUES) for (i=0; i<red->N; i++) buffer[i] = gv[i] + lv[i];
1218ac4e037SJed Brown #if !defined(PETSC_USE_COMPLEX)
1228ac4e037SJed Brown       if (imode == MAX_VALUES) for (i=0; i<red->N; i++) buffer[i] = PetscMax(gv[i],lv[i]);
1238ac4e037SJed Brown #endif
1248865f1eaSKarl Rupp     } else source = (void*)lv;
1258fa295daSBarry Smith     ierr = MPI_Reduce(source,gv,red->N,MPIU_SCALAR,(imode == ADD_VALUES) ? MPIU_SUM : MPIU_MAX,red->rank,PetscObjectComm((PetscObject)dm));CHKERRQ(ierr);
1268ac4e037SJed Brown #if !defined(PETSC_HAVE_MPI_IN_PLACE)
12799839c1eSBarry Smith     if (rank == red->rank) {ierr = PetscFree(buffer);CHKERRQ(ierr);}
1288ac4e037SJed Brown #endif
1298ac4e037SJed Brown   } break;
1308ac4e037SJed Brown   case INSERT_VALUES:
131580bdb30SBarry Smith     ierr = PetscArraycpy(gv,lv,red->n);CHKERRQ(ierr);
1328ac4e037SJed Brown     break;
133ce94432eSBarry Smith   default: SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"InsertMode not supported");
1348ac4e037SJed Brown   }
1358ac4e037SJed Brown   ierr = VecRestoreArrayRead(l,&lv);CHKERRQ(ierr);
1368ac4e037SJed Brown   ierr = VecRestoreArray(g,&gv);CHKERRQ(ierr);
1378ac4e037SJed Brown   PetscFunctionReturn(0);
1388ac4e037SJed Brown }
1398ac4e037SJed Brown 
1408ac4e037SJed Brown static PetscErrorCode DMLocalToGlobalEnd_Redundant(DM dm,Vec l,InsertMode imode,Vec g)
1418ac4e037SJed Brown {
1428ac4e037SJed Brown   PetscFunctionBegin;
1438ac4e037SJed Brown   PetscFunctionReturn(0);
1448ac4e037SJed Brown }
1458ac4e037SJed Brown 
1468ac4e037SJed Brown static PetscErrorCode DMGlobalToLocalBegin_Redundant(DM dm,Vec g,InsertMode imode,Vec l)
1478ac4e037SJed Brown {
1488ac4e037SJed Brown   PetscErrorCode    ierr;
1498ac4e037SJed Brown   DM_Redundant      *red = (DM_Redundant*)dm->data;
1508ac4e037SJed Brown   const PetscScalar *gv;
1518ac4e037SJed Brown   PetscScalar       *lv;
1528ac4e037SJed Brown 
1538ac4e037SJed Brown   PetscFunctionBegin;
1548ac4e037SJed Brown   ierr = VecGetArrayRead(g,&gv);CHKERRQ(ierr);
1558ac4e037SJed Brown   ierr = VecGetArray(l,&lv);CHKERRQ(ierr);
1568ac4e037SJed Brown   switch (imode) {
1578ac4e037SJed Brown   case INSERT_VALUES:
158580bdb30SBarry Smith     if (red->n) {ierr = PetscArraycpy(lv,gv,red->n);CHKERRQ(ierr);}
159ce94432eSBarry Smith     ierr = MPI_Bcast(lv,red->N,MPIU_SCALAR,red->rank,PetscObjectComm((PetscObject)dm));CHKERRQ(ierr);
1608ac4e037SJed Brown     break;
161ce94432eSBarry Smith   default: SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"InsertMode not supported");
1628ac4e037SJed Brown   }
1638ac4e037SJed Brown   ierr = VecRestoreArrayRead(g,&gv);CHKERRQ(ierr);
1648ac4e037SJed Brown   ierr = VecRestoreArray(l,&lv);CHKERRQ(ierr);
1658ac4e037SJed Brown   PetscFunctionReturn(0);
1668ac4e037SJed Brown }
1678ac4e037SJed Brown 
1688ac4e037SJed Brown static PetscErrorCode DMGlobalToLocalEnd_Redundant(DM dm,Vec g,InsertMode imode,Vec l)
1698ac4e037SJed Brown {
1708ac4e037SJed Brown   PetscFunctionBegin;
1718ac4e037SJed Brown   PetscFunctionReturn(0);
1728ac4e037SJed Brown }
1738ac4e037SJed Brown 
1748ac4e037SJed Brown static PetscErrorCode DMSetUp_Redundant(DM dm)
1758ac4e037SJed Brown {
1768ac4e037SJed Brown   PetscFunctionBegin;
1778ac4e037SJed Brown   PetscFunctionReturn(0);
1788ac4e037SJed Brown }
1798ac4e037SJed Brown 
1808ac4e037SJed Brown static PetscErrorCode DMView_Redundant(DM dm,PetscViewer viewer)
1818ac4e037SJed Brown {
1828ac4e037SJed Brown   PetscErrorCode ierr;
1838ac4e037SJed Brown   DM_Redundant   *red = (DM_Redundant*)dm->data;
1848ac4e037SJed Brown   PetscBool      iascii;
1858ac4e037SJed Brown 
1868ac4e037SJed Brown   PetscFunctionBegin;
187251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
1888ac4e037SJed Brown   if (iascii) {
1898ac4e037SJed Brown     ierr = PetscViewerASCIIPrintf(viewer,"redundant: rank=%D N=%D\n",red->rank,red->N);CHKERRQ(ierr);
1908ac4e037SJed Brown   }
1918ac4e037SJed Brown   PetscFunctionReturn(0);
1928ac4e037SJed Brown }
1938ac4e037SJed Brown 
194b412c318SBarry Smith static PetscErrorCode DMCreateColoring_Redundant(DM dm,ISColoringType ctype,ISColoring *coloring)
19569527181SJed Brown {
19669527181SJed Brown   DM_Redundant    *red = (DM_Redundant*)dm->data;
19769527181SJed Brown   PetscErrorCode  ierr;
19869527181SJed Brown   PetscInt        i,nloc;
19969527181SJed Brown   ISColoringValue *colors;
20069527181SJed Brown 
20169527181SJed Brown   PetscFunctionBegin;
20269527181SJed Brown   switch (ctype) {
20369527181SJed Brown   case IS_COLORING_GLOBAL:
20469527181SJed Brown     nloc = red->n;
20569527181SJed Brown     break;
2065bdb020cSBarry Smith   case IS_COLORING_LOCAL:
20769527181SJed Brown     nloc = red->N;
20869527181SJed Brown     break;
209ce94432eSBarry Smith   default: SETERRQ1(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONG,"Unknown ISColoringType %d",(int)ctype);
21069527181SJed Brown   }
211785e854fSJed Brown   ierr = PetscMalloc1(nloc,&colors);CHKERRQ(ierr);
21269527181SJed Brown   for (i=0; i<nloc; i++) colors[i] = i;
213aaf3ff59SMatthew G. Knepley   ierr = ISColoringCreate(PetscObjectComm((PetscObject)dm),red->N,nloc,colors,PETSC_OWN_POINTER,coloring);CHKERRQ(ierr);
21469527181SJed Brown   ierr = ISColoringSetType(*coloring,ctype);CHKERRQ(ierr);
21569527181SJed Brown   PetscFunctionReturn(0);
21669527181SJed Brown }
2178ac4e037SJed Brown 
2188ac4e037SJed Brown static PetscErrorCode DMRefine_Redundant(DM dmc,MPI_Comm comm,DM *dmf)
2198ac4e037SJed Brown {
2208ac4e037SJed Brown   PetscErrorCode ierr;
2218ac4e037SJed Brown   PetscMPIInt    flag;
2228ac4e037SJed Brown   DM_Redundant   *redc = (DM_Redundant*)dmc->data;
2238ac4e037SJed Brown 
2248ac4e037SJed Brown   PetscFunctionBegin;
225ce94432eSBarry Smith   if (comm == MPI_COMM_NULL) {
226ce94432eSBarry Smith     ierr = PetscObjectGetComm((PetscObject)dmc,&comm);CHKERRQ(ierr);
227ce94432eSBarry Smith   }
228ce94432eSBarry Smith   ierr = MPI_Comm_compare(PetscObjectComm((PetscObject)dmc),comm,&flag);CHKERRQ(ierr);
229ce94432eSBarry Smith   if (flag != MPI_CONGRUENT && flag != MPI_IDENT) SETERRQ(PetscObjectComm((PetscObject)dmc),PETSC_ERR_SUP,"cannot change communicators");
2308ac4e037SJed Brown   ierr = DMRedundantCreate(comm,redc->rank,redc->N,dmf);CHKERRQ(ierr);
2318ac4e037SJed Brown   PetscFunctionReturn(0);
2328ac4e037SJed Brown }
2338ac4e037SJed Brown 
2348ac4e037SJed Brown static PetscErrorCode DMCoarsen_Redundant(DM dmf,MPI_Comm comm,DM *dmc)
2358ac4e037SJed Brown {
2368ac4e037SJed Brown   PetscErrorCode ierr;
2378ac4e037SJed Brown   PetscMPIInt    flag;
2388ac4e037SJed Brown   DM_Redundant   *redf = (DM_Redundant*)dmf->data;
2398ac4e037SJed Brown 
2408ac4e037SJed Brown   PetscFunctionBegin;
241ce94432eSBarry Smith   if (comm == MPI_COMM_NULL) {
242ce94432eSBarry Smith     ierr = PetscObjectGetComm((PetscObject)dmf,&comm);CHKERRQ(ierr);
243ce94432eSBarry Smith   }
244ce94432eSBarry Smith   ierr = MPI_Comm_compare(PetscObjectComm((PetscObject)dmf),comm,&flag);CHKERRQ(ierr);
245ce94432eSBarry Smith   if (flag != MPI_CONGRUENT && flag != MPI_IDENT) SETERRQ(PetscObjectComm((PetscObject)dmf),PETSC_ERR_SUP,"cannot change communicators");
2468ac4e037SJed Brown   ierr = DMRedundantCreate(comm,redf->rank,redf->N,dmc);CHKERRQ(ierr);
2478ac4e037SJed Brown   PetscFunctionReturn(0);
2488ac4e037SJed Brown }
2498ac4e037SJed Brown 
250e727c939SJed Brown static PetscErrorCode DMCreateInterpolation_Redundant(DM dmc,DM dmf,Mat *P,Vec *scale)
2518ac4e037SJed Brown {
2528ac4e037SJed Brown   PetscErrorCode ierr;
2538ac4e037SJed Brown   DM_Redundant   *redc = (DM_Redundant*)dmc->data;
2548ac4e037SJed Brown   DM_Redundant   *redf = (DM_Redundant*)dmf->data;
2558ac4e037SJed Brown   PetscMPIInt    flag;
2568ac4e037SJed Brown   PetscInt       i,rstart,rend;
2578ac4e037SJed Brown 
2588ac4e037SJed Brown   PetscFunctionBegin;
259ce94432eSBarry Smith   ierr = MPI_Comm_compare(PetscObjectComm((PetscObject)dmc),PetscObjectComm((PetscObject)dmf),&flag);CHKERRQ(ierr);
260ce94432eSBarry Smith   if (flag != MPI_CONGRUENT && flag != MPI_IDENT) SETERRQ(PetscObjectComm((PetscObject)dmf),PETSC_ERR_SUP,"cannot change communicators");
261ce94432eSBarry Smith   if (redc->rank != redf->rank) SETERRQ(PetscObjectComm((PetscObject)dmf),PETSC_ERR_ARG_INCOMP,"Owning rank does not match");
262ce94432eSBarry Smith   if (redc->N != redf->N) SETERRQ(PetscObjectComm((PetscObject)dmf),PETSC_ERR_ARG_INCOMP,"Global size does not match");
263ce94432eSBarry Smith   ierr = MatCreate(PetscObjectComm((PetscObject)dmc),P);CHKERRQ(ierr);
2648ac4e037SJed Brown   ierr = MatSetSizes(*P,redc->n,redc->n,redc->N,redc->N);CHKERRQ(ierr);
2658ac4e037SJed Brown   ierr = MatSetType(*P,MATAIJ);CHKERRQ(ierr);
2668ac4e037SJed Brown   ierr = MatSeqAIJSetPreallocation(*P,1,0);CHKERRQ(ierr);
2678ac4e037SJed Brown   ierr = MatMPIAIJSetPreallocation(*P,1,0,0,0);CHKERRQ(ierr);
2688ac4e037SJed Brown   ierr = MatGetOwnershipRange(*P,&rstart,&rend);CHKERRQ(ierr);
2698ac4e037SJed Brown   for (i=rstart; i<rend; i++) {ierr = MatSetValue(*P,i,i,1.0,INSERT_VALUES);CHKERRQ(ierr);}
2708ac4e037SJed Brown   ierr = MatAssemblyBegin(*P,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
2718ac4e037SJed Brown   ierr = MatAssemblyEnd(*P,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
272e727c939SJed Brown   if (scale) {ierr = DMCreateInterpolationScale(dmc,dmf,*P,scale);CHKERRQ(ierr);}
2738ac4e037SJed Brown   PetscFunctionReturn(0);
2748ac4e037SJed Brown }
2758ac4e037SJed Brown 
2768ac4e037SJed Brown /*@
2778ac4e037SJed Brown     DMRedundantSetSize - Sets the size of a densely coupled redundant object
2788ac4e037SJed Brown 
279d083f849SBarry Smith     Collective on dm
2808ac4e037SJed Brown 
2818ac4e037SJed Brown     Input Parameter:
2828ac4e037SJed Brown +   dm - redundant DM
2838ac4e037SJed Brown .   rank - rank of process to own redundant degrees of freedom
2848ac4e037SJed Brown -   N - total number of redundant degrees of freedom
2858ac4e037SJed Brown 
2868ac4e037SJed Brown     Level: advanced
2878ac4e037SJed Brown 
2888ac4e037SJed Brown .seealso DMDestroy(), DMCreateGlobalVector(), DMRedundantCreate(), DMRedundantGetSize()
2898ac4e037SJed Brown @*/
290907376e6SBarry Smith PetscErrorCode DMRedundantSetSize(DM dm,PetscMPIInt rank,PetscInt N)
2918ac4e037SJed Brown {
2928ac4e037SJed Brown   PetscErrorCode ierr;
2938ac4e037SJed Brown 
2948ac4e037SJed Brown   PetscFunctionBegin;
2958ac4e037SJed Brown   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
2968ac4e037SJed Brown   PetscValidType(dm,1);
29769fbec6eSBarry Smith   PetscValidLogicalCollectiveMPIInt(dm,rank,2);
2988ac4e037SJed Brown   PetscValidLogicalCollectiveInt(dm,N,3);
29983515869SBarry Smith   ierr = PetscTryMethod(dm,"DMRedundantSetSize_C",(DM,PetscMPIInt,PetscInt),(dm,rank,N));CHKERRQ(ierr);
3008ac4e037SJed Brown   PetscFunctionReturn(0);
3018ac4e037SJed Brown }
3028ac4e037SJed Brown 
3038ac4e037SJed Brown /*@
3048ac4e037SJed Brown     DMRedundantGetSize - Gets the size of a densely coupled redundant object
3058ac4e037SJed Brown 
3068ac4e037SJed Brown     Not Collective
3078ac4e037SJed Brown 
3088ac4e037SJed Brown     Input Parameter:
309a2b725a8SWilliam Gropp .   dm - redundant DM
3108ac4e037SJed Brown 
3118ac4e037SJed Brown     Output Parameters:
3120298fd71SBarry Smith +   rank - rank of process to own redundant degrees of freedom (or NULL)
3130298fd71SBarry Smith -   N - total number of redundant degrees of freedom (or NULL)
3148ac4e037SJed Brown 
3158ac4e037SJed Brown     Level: advanced
3168ac4e037SJed Brown 
3178ac4e037SJed Brown .seealso DMDestroy(), DMCreateGlobalVector(), DMRedundantCreate(), DMRedundantSetSize()
3188ac4e037SJed Brown @*/
319907376e6SBarry Smith PetscErrorCode DMRedundantGetSize(DM dm,PetscMPIInt *rank,PetscInt *N)
3208ac4e037SJed Brown {
3218ac4e037SJed Brown   PetscErrorCode ierr;
3228ac4e037SJed Brown 
3238ac4e037SJed Brown   PetscFunctionBegin;
3248ac4e037SJed Brown   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
3258ac4e037SJed Brown   PetscValidType(dm,1);
32683515869SBarry Smith   ierr = PetscUseMethod(dm,"DMRedundantGetSize_C",(DM,PetscMPIInt*,PetscInt*),(dm,rank,N));CHKERRQ(ierr);
3278ac4e037SJed Brown   PetscFunctionReturn(0);
3288ac4e037SJed Brown }
3298ac4e037SJed Brown 
330907376e6SBarry Smith static PetscErrorCode DMRedundantSetSize_Redundant(DM dm,PetscMPIInt rank,PetscInt N)
3318ac4e037SJed Brown {
3328ac4e037SJed Brown   DM_Redundant   *red = (DM_Redundant*)dm->data;
3338ac4e037SJed Brown   PetscErrorCode ierr;
3348ac4e037SJed Brown   PetscMPIInt    myrank;
335*73c9e547SStefano Zampini   PetscInt       i,*globals;
3368ac4e037SJed Brown 
3378ac4e037SJed Brown   PetscFunctionBegin;
338ce94432eSBarry Smith   ierr      = MPI_Comm_rank(PetscObjectComm((PetscObject)dm),&myrank);CHKERRQ(ierr);
3398ac4e037SJed Brown   red->rank = rank;
3408ac4e037SJed Brown   red->N    = N;
3418ac4e037SJed Brown   red->n    = (myrank == rank) ? N : 0;
342*73c9e547SStefano Zampini 
343*73c9e547SStefano Zampini   /* mapping is setup here */
344*73c9e547SStefano Zampini   ierr = PetscMalloc1(red->N,&globals);CHKERRQ(ierr);
345*73c9e547SStefano Zampini   for (i=0; i<red->N; i++) globals[i] = i;
346*73c9e547SStefano Zampini   ierr = ISLocalToGlobalMappingDestroy(&dm->ltogmap);CHKERRQ(ierr);
347*73c9e547SStefano Zampini   ierr = ISLocalToGlobalMappingCreate(PetscObjectComm((PetscObject)dm),1,red->N,globals,PETSC_OWN_POINTER,&dm->ltogmap);CHKERRQ(ierr);
3488ac4e037SJed Brown   PetscFunctionReturn(0);
3498ac4e037SJed Brown }
3508ac4e037SJed Brown 
3511e6b0712SBarry Smith static PetscErrorCode DMRedundantGetSize_Redundant(DM dm,PetscInt *rank,PetscInt *N)
3528ac4e037SJed Brown {
3538ac4e037SJed Brown   DM_Redundant *red = (DM_Redundant*)dm->data;
3548ac4e037SJed Brown 
3558ac4e037SJed Brown   PetscFunctionBegin;
3568ac4e037SJed Brown   if (rank) *rank = red->rank;
3578ac4e037SJed Brown   if (N)    *N = red->N;
3588ac4e037SJed Brown   PetscFunctionReturn(0);
3598ac4e037SJed Brown }
3608ac4e037SJed Brown 
36136623e74SStefano Zampini static PetscErrorCode DMSetUpGLVisViewer_Redundant(PetscObject odm, PetscViewer viewer)
36236623e74SStefano Zampini {
36336623e74SStefano Zampini   PetscFunctionBegin;
36436623e74SStefano Zampini   PetscFunctionReturn(0);
36536623e74SStefano Zampini }
36636623e74SStefano Zampini 
3673efe6655SBarry Smith /*MC
3683efe6655SBarry Smith    DMREDUNDANT = "redundant" - A DM object that is used to manage data for a small set of dense globally coupled variables.
3693efe6655SBarry Smith          In the global representation of the vector the variables are all stored on a single MPI process (all the other MPI processes have
3703efe6655SBarry Smith          no variables) in the local representation all the variables are stored on ALL the MPI processes (because they are all needed for each
3713efe6655SBarry Smith          processes local computations).
3723efe6655SBarry Smith 
3733efe6655SBarry Smith          This DM is generally used inside a DMCOMPOSITE object. For example, it may be used to store continuation parameters for a bifurcation problem.
3743efe6655SBarry Smith 
3753efe6655SBarry Smith 
3763efe6655SBarry Smith   Level: intermediate
3773efe6655SBarry Smith 
3781abcec8cSBarry Smith .seealso: DMType, DMCOMPOSITE,  DMCreate(), DMRedundantSetSize(), DMRedundantGetSize()
3793efe6655SBarry Smith M*/
3803efe6655SBarry Smith 
3818cc058d9SJed Brown PETSC_EXTERN PetscErrorCode DMCreate_Redundant(DM dm)
3828ac4e037SJed Brown {
3838ac4e037SJed Brown   PetscErrorCode ierr;
3848ac4e037SJed Brown   DM_Redundant   *red;
3858ac4e037SJed Brown 
3868ac4e037SJed Brown   PetscFunctionBegin;
387b00a9115SJed Brown   ierr     = PetscNewLog(dm,&red);CHKERRQ(ierr);
3888ac4e037SJed Brown   dm->data = red;
3898ac4e037SJed Brown 
3908ac4e037SJed Brown   dm->ops->setup               = DMSetUp_Redundant;
3918ac4e037SJed Brown   dm->ops->view                = DMView_Redundant;
3928ac4e037SJed Brown   dm->ops->createglobalvector  = DMCreateGlobalVector_Redundant;
3938ac4e037SJed Brown   dm->ops->createlocalvector   = DMCreateLocalVector_Redundant;
39425296bd5SBarry Smith   dm->ops->creatematrix        = DMCreateMatrix_Redundant;
3958ac4e037SJed Brown   dm->ops->destroy             = DMDestroy_Redundant;
3968ac4e037SJed Brown   dm->ops->globaltolocalbegin  = DMGlobalToLocalBegin_Redundant;
3978ac4e037SJed Brown   dm->ops->globaltolocalend    = DMGlobalToLocalEnd_Redundant;
3988ac4e037SJed Brown   dm->ops->localtoglobalbegin  = DMLocalToGlobalBegin_Redundant;
3998ac4e037SJed Brown   dm->ops->localtoglobalend    = DMLocalToGlobalEnd_Redundant;
4008ac4e037SJed Brown   dm->ops->refine              = DMRefine_Redundant;
4018ac4e037SJed Brown   dm->ops->coarsen             = DMCoarsen_Redundant;
40225296bd5SBarry Smith   dm->ops->createinterpolation = DMCreateInterpolation_Redundant;
403e727c939SJed Brown   dm->ops->getcoloring         = DMCreateColoring_Redundant;
4048865f1eaSKarl Rupp 
405bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)dm,"DMRedundantSetSize_C",DMRedundantSetSize_Redundant);CHKERRQ(ierr);
406bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)dm,"DMRedundantGetSize_C",DMRedundantGetSize_Redundant);CHKERRQ(ierr);
40736623e74SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)dm,"DMSetUpGLVisViewer_C",DMSetUpGLVisViewer_Redundant);CHKERRQ(ierr);
4088ac4e037SJed Brown   PetscFunctionReturn(0);
4098ac4e037SJed Brown }
4108ac4e037SJed Brown 
4118ac4e037SJed Brown /*@C
4128ac4e037SJed Brown     DMRedundantCreate - Creates a DM object, used to manage data for dense globally coupled variables
4138ac4e037SJed Brown 
414d083f849SBarry Smith     Collective
4158ac4e037SJed Brown 
4168ac4e037SJed Brown     Input Parameter:
4178ac4e037SJed Brown +   comm - the processors that will share the global vector
4188ac4e037SJed Brown .   rank - rank to own the redundant values
4198ac4e037SJed Brown -   N - total number of degrees of freedom
4208ac4e037SJed Brown 
4218ac4e037SJed Brown     Output Parameters:
422907376e6SBarry Smith .   dm - the redundant DM
4238ac4e037SJed Brown 
4248ac4e037SJed Brown     Level: advanced
4258ac4e037SJed Brown 
4263efe6655SBarry Smith .seealso DMDestroy(), DMCreateGlobalVector(), DMCreateMatrix(), DMCompositeAddDM(), DMREDUNDANT, DMSetType(), DMRedundantSetSize(), DMRedundantGetSize()
4278ac4e037SJed Brown 
4288ac4e037SJed Brown @*/
429907376e6SBarry Smith PetscErrorCode DMRedundantCreate(MPI_Comm comm,PetscMPIInt rank,PetscInt N,DM *dm)
4308ac4e037SJed Brown {
4318ac4e037SJed Brown   PetscErrorCode ierr;
4328ac4e037SJed Brown 
4338ac4e037SJed Brown   PetscFunctionBegin;
4348ac4e037SJed Brown   PetscValidPointer(dm,2);
4358ac4e037SJed Brown   ierr = DMCreate(comm,dm);CHKERRQ(ierr);
4368ac4e037SJed Brown   ierr = DMSetType(*dm,DMREDUNDANT);CHKERRQ(ierr);
4378ac4e037SJed Brown   ierr = DMRedundantSetSize(*dm,rank,N);CHKERRQ(ierr);
4388ac4e037SJed Brown   ierr = DMSetUp(*dm);CHKERRQ(ierr);
4398ac4e037SJed Brown   PetscFunctionReturn(0);
4408ac4e037SJed Brown }
441