xref: /petsc/src/dm/impls/redundant/dmredundant.c (revision 1dca8a0504492127e77eac64bc165d7372dd6d63)
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;
1345b6f7e9SBarry Smith   ISLocalToGlobalMapping ltog;
14d2e2b695SJed Brown   PetscInt               i,rstart,rend,*cols;
15d2e2b695SJed Brown   PetscScalar            *vals;
168ac4e037SJed Brown 
178ac4e037SJed Brown   PetscFunctionBegin;
189566063dSJacob Faibussowitsch   PetscCall(MatCreate(PetscObjectComm((PetscObject)dm),J));
199566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(*J,red->n,red->n,red->N,red->N));
209566063dSJacob Faibussowitsch   PetscCall(MatSetType(*J,dm->mattype));
219566063dSJacob Faibussowitsch   PetscCall(MatSeqAIJSetPreallocation(*J,red->n,NULL));
229566063dSJacob Faibussowitsch   PetscCall(MatSeqBAIJSetPreallocation(*J,1,red->n,NULL));
239566063dSJacob Faibussowitsch   PetscCall(MatMPIAIJSetPreallocation(*J,red->n,NULL,red->N-red->n,NULL));
249566063dSJacob Faibussowitsch   PetscCall(MatMPIBAIJSetPreallocation(*J,1,red->n,NULL,red->N-red->n,NULL));
258ac4e037SJed Brown 
269566063dSJacob Faibussowitsch   PetscCall(DMGetLocalToGlobalMapping(dm,&ltog));
279566063dSJacob Faibussowitsch   PetscCall(MatSetLocalToGlobalMapping(*J,ltog,ltog));
289566063dSJacob Faibussowitsch   PetscCall(MatSetDM(*J,dm));
29d2e2b695SJed Brown 
309566063dSJacob Faibussowitsch   PetscCall(PetscMalloc2(red->N,&cols,red->N,&vals));
31d2e2b695SJed Brown   for (i=0; i<red->N; i++) {
32d2e2b695SJed Brown     cols[i] = i;
33d2e2b695SJed Brown     vals[i] = 0.0;
34d2e2b695SJed Brown   }
359566063dSJacob Faibussowitsch   PetscCall(MatGetOwnershipRange(*J,&rstart,&rend));
36d2e2b695SJed Brown   for (i=rstart; i<rend; i++) {
379566063dSJacob Faibussowitsch     PetscCall(MatSetValues(*J,1,&i,red->N,cols,vals,INSERT_VALUES));
38d2e2b695SJed Brown   }
399566063dSJacob Faibussowitsch   PetscCall(PetscFree2(cols,vals));
409566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(*J,MAT_FINAL_ASSEMBLY));
419566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(*J,MAT_FINAL_ASSEMBLY));
428ac4e037SJed Brown   PetscFunctionReturn(0);
438ac4e037SJed Brown }
448ac4e037SJed Brown 
458ac4e037SJed Brown static PetscErrorCode DMDestroy_Redundant(DM dm)
468ac4e037SJed Brown {
478ac4e037SJed Brown   PetscFunctionBegin;
489566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)dm,"DMRedundantSetSize_C",NULL));
499566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)dm,"DMRedundantGetSize_C",NULL));
509566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)dm,"DMSetUpGLVisViewer_C",NULL));
51435a35e8SMatthew G Knepley   /* This was originally freed in DMDestroy(), but that prevents reference counting of backend objects */
529566063dSJacob Faibussowitsch   PetscCall(PetscFree(dm->data));
538ac4e037SJed Brown   PetscFunctionReturn(0);
548ac4e037SJed Brown }
558ac4e037SJed Brown 
568ac4e037SJed Brown static PetscErrorCode DMCreateGlobalVector_Redundant(DM dm,Vec *gvec)
578ac4e037SJed Brown {
588ac4e037SJed Brown   DM_Redundant           *red = (DM_Redundant*)dm->data;
5945b6f7e9SBarry Smith   ISLocalToGlobalMapping ltog;
608ac4e037SJed Brown 
618ac4e037SJed Brown   PetscFunctionBegin;
628ac4e037SJed Brown   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
638ac4e037SJed Brown   PetscValidPointer(gvec,2);
64ea78f98cSLisandro Dalcin   *gvec = NULL;
659566063dSJacob Faibussowitsch   PetscCall(VecCreate(PetscObjectComm((PetscObject)dm),gvec));
669566063dSJacob Faibussowitsch   PetscCall(VecSetSizes(*gvec,red->n,red->N));
679566063dSJacob Faibussowitsch   PetscCall(VecSetType(*gvec,dm->vectype));
689566063dSJacob Faibussowitsch   PetscCall(DMGetLocalToGlobalMapping(dm,&ltog));
699566063dSJacob Faibussowitsch   PetscCall(VecSetLocalToGlobalMapping(*gvec,ltog));
709566063dSJacob Faibussowitsch   PetscCall(VecSetDM(*gvec,dm));
718ac4e037SJed Brown   PetscFunctionReturn(0);
728ac4e037SJed Brown }
738ac4e037SJed Brown 
748ac4e037SJed Brown static PetscErrorCode DMCreateLocalVector_Redundant(DM dm,Vec *lvec)
758ac4e037SJed Brown {
768ac4e037SJed Brown   DM_Redundant   *red = (DM_Redundant*)dm->data;
778ac4e037SJed Brown 
788ac4e037SJed Brown   PetscFunctionBegin;
798ac4e037SJed Brown   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
808ac4e037SJed Brown   PetscValidPointer(lvec,2);
81ea78f98cSLisandro Dalcin   *lvec = NULL;
829566063dSJacob Faibussowitsch   PetscCall(VecCreate(PETSC_COMM_SELF,lvec));
839566063dSJacob Faibussowitsch   PetscCall(VecSetSizes(*lvec,red->N,red->N));
849566063dSJacob Faibussowitsch   PetscCall(VecSetType(*lvec,dm->vectype));
859566063dSJacob Faibussowitsch   PetscCall(VecSetDM(*lvec,dm));
868ac4e037SJed Brown   PetscFunctionReturn(0);
878ac4e037SJed Brown }
888ac4e037SJed Brown 
898ac4e037SJed Brown static PetscErrorCode DMLocalToGlobalBegin_Redundant(DM dm,Vec l,InsertMode imode,Vec g)
908ac4e037SJed Brown {
918ac4e037SJed Brown   DM_Redundant      *red = (DM_Redundant*)dm->data;
928ac4e037SJed Brown   const PetscScalar *lv;
938ac4e037SJed Brown   PetscScalar       *gv;
948ac4e037SJed Brown   PetscMPIInt       rank;
958ac4e037SJed Brown 
968ac4e037SJed Brown   PetscFunctionBegin;
979566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)dm),&rank));
989566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(l,&lv));
999566063dSJacob Faibussowitsch   PetscCall(VecGetArray(g,&gv));
1008ac4e037SJed Brown   switch (imode) {
1018ac4e037SJed Brown   case ADD_VALUES:
1028ac4e037SJed Brown   case MAX_VALUES:
1038ac4e037SJed Brown   {
1048ac4e037SJed Brown     void        *source;
1058ac4e037SJed Brown     PetscScalar *buffer;
1068ac4e037SJed Brown     PetscInt    i;
1078ac4e037SJed Brown     if (rank == red->rank) {
1088ac4e037SJed Brown       buffer = gv;
1098ac4e037SJed Brown       source = MPI_IN_PLACE;
1108ac4e037SJed Brown       if (imode == ADD_VALUES) for (i=0; i<red->N; i++) buffer[i] = gv[i] + lv[i];
1118ac4e037SJed Brown #if !defined(PETSC_USE_COMPLEX)
1128ac4e037SJed Brown       if (imode == MAX_VALUES) for (i=0; i<red->N; i++) buffer[i] = PetscMax(gv[i],lv[i]);
1138ac4e037SJed Brown #endif
1148865f1eaSKarl Rupp     } else source = (void*)lv;
1159566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Reduce(source,gv,red->N,MPIU_SCALAR,(imode == ADD_VALUES) ? MPIU_SUM : MPIU_MAX,red->rank,PetscObjectComm((PetscObject)dm)));
1168ac4e037SJed Brown   } break;
1178ac4e037SJed Brown   case INSERT_VALUES:
1189566063dSJacob Faibussowitsch     PetscCall(PetscArraycpy(gv,lv,red->n));
1198ac4e037SJed Brown     break;
120ce94432eSBarry Smith   default: SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"InsertMode not supported");
1218ac4e037SJed Brown   }
1229566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(l,&lv));
1239566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(g,&gv));
1248ac4e037SJed Brown   PetscFunctionReturn(0);
1258ac4e037SJed Brown }
1268ac4e037SJed Brown 
1278ac4e037SJed Brown static PetscErrorCode DMLocalToGlobalEnd_Redundant(DM dm,Vec l,InsertMode imode,Vec g)
1288ac4e037SJed Brown {
1298ac4e037SJed Brown   PetscFunctionBegin;
1308ac4e037SJed Brown   PetscFunctionReturn(0);
1318ac4e037SJed Brown }
1328ac4e037SJed Brown 
1338ac4e037SJed Brown static PetscErrorCode DMGlobalToLocalBegin_Redundant(DM dm,Vec g,InsertMode imode,Vec l)
1348ac4e037SJed Brown {
1358ac4e037SJed Brown   DM_Redundant      *red = (DM_Redundant*)dm->data;
1368ac4e037SJed Brown   const PetscScalar *gv;
1378ac4e037SJed Brown   PetscScalar       *lv;
1388ac4e037SJed Brown 
1398ac4e037SJed Brown   PetscFunctionBegin;
1409566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(g,&gv));
1419566063dSJacob Faibussowitsch   PetscCall(VecGetArray(l,&lv));
1428ac4e037SJed Brown   switch (imode) {
1438ac4e037SJed Brown   case INSERT_VALUES:
1449566063dSJacob Faibussowitsch     if (red->n) PetscCall(PetscArraycpy(lv,gv,red->n));
1459566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Bcast(lv,red->N,MPIU_SCALAR,red->rank,PetscObjectComm((PetscObject)dm)));
1468ac4e037SJed Brown     break;
147ce94432eSBarry Smith   default: SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"InsertMode not supported");
1488ac4e037SJed Brown   }
1499566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(g,&gv));
1509566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(l,&lv));
1518ac4e037SJed Brown   PetscFunctionReturn(0);
1528ac4e037SJed Brown }
1538ac4e037SJed Brown 
1548ac4e037SJed Brown static PetscErrorCode DMGlobalToLocalEnd_Redundant(DM dm,Vec g,InsertMode imode,Vec l)
1558ac4e037SJed Brown {
1568ac4e037SJed Brown   PetscFunctionBegin;
1578ac4e037SJed Brown   PetscFunctionReturn(0);
1588ac4e037SJed Brown }
1598ac4e037SJed Brown 
1608ac4e037SJed Brown static PetscErrorCode DMSetUp_Redundant(DM dm)
1618ac4e037SJed Brown {
1628ac4e037SJed Brown   PetscFunctionBegin;
1638ac4e037SJed Brown   PetscFunctionReturn(0);
1648ac4e037SJed Brown }
1658ac4e037SJed Brown 
1668ac4e037SJed Brown static PetscErrorCode DMView_Redundant(DM dm,PetscViewer viewer)
1678ac4e037SJed Brown {
1688ac4e037SJed Brown   DM_Redundant   *red = (DM_Redundant*)dm->data;
1698ac4e037SJed Brown   PetscBool      iascii;
1708ac4e037SJed Brown 
1718ac4e037SJed Brown   PetscFunctionBegin;
1729566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii));
1738ac4e037SJed Brown   if (iascii) {
17463a3b9bcSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer,"redundant: rank=%d N=%" PetscInt_FMT "\n",red->rank,red->N));
1758ac4e037SJed Brown   }
1768ac4e037SJed Brown   PetscFunctionReturn(0);
1778ac4e037SJed Brown }
1788ac4e037SJed Brown 
179b412c318SBarry Smith static PetscErrorCode DMCreateColoring_Redundant(DM dm,ISColoringType ctype,ISColoring *coloring)
18069527181SJed Brown {
18169527181SJed Brown   DM_Redundant    *red = (DM_Redundant*)dm->data;
18269527181SJed Brown   PetscInt        i,nloc;
18369527181SJed Brown   ISColoringValue *colors;
18469527181SJed Brown 
18569527181SJed Brown   PetscFunctionBegin;
18669527181SJed Brown   switch (ctype) {
18769527181SJed Brown   case IS_COLORING_GLOBAL:
18869527181SJed Brown     nloc = red->n;
18969527181SJed Brown     break;
1905bdb020cSBarry Smith   case IS_COLORING_LOCAL:
19169527181SJed Brown     nloc = red->N;
19269527181SJed Brown     break;
19398921bdaSJacob Faibussowitsch   default: SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONG,"Unknown ISColoringType %d",(int)ctype);
19469527181SJed Brown   }
1959566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(nloc,&colors));
19669527181SJed Brown   for (i=0; i<nloc; i++) colors[i] = i;
1979566063dSJacob Faibussowitsch   PetscCall(ISColoringCreate(PetscObjectComm((PetscObject)dm),red->N,nloc,colors,PETSC_OWN_POINTER,coloring));
1989566063dSJacob Faibussowitsch   PetscCall(ISColoringSetType(*coloring,ctype));
19969527181SJed Brown   PetscFunctionReturn(0);
20069527181SJed Brown }
2018ac4e037SJed Brown 
2028ac4e037SJed Brown static PetscErrorCode DMRefine_Redundant(DM dmc,MPI_Comm comm,DM *dmf)
2038ac4e037SJed Brown {
2048ac4e037SJed Brown   PetscMPIInt    flag;
2058ac4e037SJed Brown   DM_Redundant   *redc = (DM_Redundant*)dmc->data;
2068ac4e037SJed Brown 
2078ac4e037SJed Brown   PetscFunctionBegin;
208ce94432eSBarry Smith   if (comm == MPI_COMM_NULL) {
2099566063dSJacob Faibussowitsch     PetscCall(PetscObjectGetComm((PetscObject)dmc,&comm));
210ce94432eSBarry Smith   }
2119566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_compare(PetscObjectComm((PetscObject)dmc),comm,&flag));
212*1dca8a05SBarry Smith   PetscCheck(flag == MPI_CONGRUENT || flag == MPI_IDENT,PetscObjectComm((PetscObject)dmc),PETSC_ERR_SUP,"cannot change communicators");
2139566063dSJacob Faibussowitsch   PetscCall(DMRedundantCreate(comm,redc->rank,redc->N,dmf));
2148ac4e037SJed Brown   PetscFunctionReturn(0);
2158ac4e037SJed Brown }
2168ac4e037SJed Brown 
2178ac4e037SJed Brown static PetscErrorCode DMCoarsen_Redundant(DM dmf,MPI_Comm comm,DM *dmc)
2188ac4e037SJed Brown {
2198ac4e037SJed Brown   PetscMPIInt    flag;
2208ac4e037SJed Brown   DM_Redundant   *redf = (DM_Redundant*)dmf->data;
2218ac4e037SJed Brown 
2228ac4e037SJed Brown   PetscFunctionBegin;
223ce94432eSBarry Smith   if (comm == MPI_COMM_NULL) {
2249566063dSJacob Faibussowitsch     PetscCall(PetscObjectGetComm((PetscObject)dmf,&comm));
225ce94432eSBarry Smith   }
2269566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_compare(PetscObjectComm((PetscObject)dmf),comm,&flag));
227*1dca8a05SBarry Smith   PetscCheck(flag == MPI_CONGRUENT || flag == MPI_IDENT,PetscObjectComm((PetscObject)dmf),PETSC_ERR_SUP,"cannot change communicators");
2289566063dSJacob Faibussowitsch   PetscCall(DMRedundantCreate(comm,redf->rank,redf->N,dmc));
2298ac4e037SJed Brown   PetscFunctionReturn(0);
2308ac4e037SJed Brown }
2318ac4e037SJed Brown 
232e727c939SJed Brown static PetscErrorCode DMCreateInterpolation_Redundant(DM dmc,DM dmf,Mat *P,Vec *scale)
2338ac4e037SJed Brown {
2348ac4e037SJed Brown   DM_Redundant   *redc = (DM_Redundant*)dmc->data;
2358ac4e037SJed Brown   DM_Redundant   *redf = (DM_Redundant*)dmf->data;
2368ac4e037SJed Brown   PetscMPIInt    flag;
2378ac4e037SJed Brown   PetscInt       i,rstart,rend;
2388ac4e037SJed Brown 
2398ac4e037SJed Brown   PetscFunctionBegin;
2409566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_compare(PetscObjectComm((PetscObject)dmc),PetscObjectComm((PetscObject)dmf),&flag));
241*1dca8a05SBarry Smith   PetscCheck(flag == MPI_CONGRUENT || flag == MPI_IDENT,PetscObjectComm((PetscObject)dmf),PETSC_ERR_SUP,"cannot change communicators");
24208401ef6SPierre Jolivet   PetscCheck(redc->rank == redf->rank,PetscObjectComm((PetscObject)dmf),PETSC_ERR_ARG_INCOMP,"Owning rank does not match");
24308401ef6SPierre Jolivet   PetscCheck(redc->N == redf->N,PetscObjectComm((PetscObject)dmf),PETSC_ERR_ARG_INCOMP,"Global size does not match");
2449566063dSJacob Faibussowitsch   PetscCall(MatCreate(PetscObjectComm((PetscObject)dmc),P));
2459566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(*P,redc->n,redc->n,redc->N,redc->N));
2469566063dSJacob Faibussowitsch   PetscCall(MatSetType(*P,MATAIJ));
2479566063dSJacob Faibussowitsch   PetscCall(MatSeqAIJSetPreallocation(*P,1,NULL));
2489566063dSJacob Faibussowitsch   PetscCall(MatMPIAIJSetPreallocation(*P,1,NULL,0,NULL));
2499566063dSJacob Faibussowitsch   PetscCall(MatGetOwnershipRange(*P,&rstart,&rend));
2509566063dSJacob Faibussowitsch   for (i=rstart; i<rend; i++) PetscCall(MatSetValue(*P,i,i,1.0,INSERT_VALUES));
2519566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(*P,MAT_FINAL_ASSEMBLY));
2529566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(*P,MAT_FINAL_ASSEMBLY));
2539566063dSJacob Faibussowitsch   if (scale) PetscCall(DMCreateInterpolationScale(dmc,dmf,*P,scale));
2548ac4e037SJed Brown   PetscFunctionReturn(0);
2558ac4e037SJed Brown }
2568ac4e037SJed Brown 
2578ac4e037SJed Brown /*@
2588ac4e037SJed Brown     DMRedundantSetSize - Sets the size of a densely coupled redundant object
2598ac4e037SJed Brown 
260d083f849SBarry Smith     Collective on dm
2618ac4e037SJed Brown 
262d8d19677SJose E. Roman     Input Parameters:
2638ac4e037SJed Brown +   dm - redundant DM
2648ac4e037SJed Brown .   rank - rank of process to own redundant degrees of freedom
2658ac4e037SJed Brown -   N - total number of redundant degrees of freedom
2668ac4e037SJed Brown 
2678ac4e037SJed Brown     Level: advanced
2688ac4e037SJed Brown 
2698ac4e037SJed Brown .seealso DMDestroy(), DMCreateGlobalVector(), DMRedundantCreate(), DMRedundantGetSize()
2708ac4e037SJed Brown @*/
271907376e6SBarry Smith PetscErrorCode DMRedundantSetSize(DM dm,PetscMPIInt rank,PetscInt N)
2728ac4e037SJed Brown {
2738ac4e037SJed Brown   PetscFunctionBegin;
2748ac4e037SJed Brown   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
2758ac4e037SJed Brown   PetscValidType(dm,1);
27669fbec6eSBarry Smith   PetscValidLogicalCollectiveMPIInt(dm,rank,2);
2778ac4e037SJed Brown   PetscValidLogicalCollectiveInt(dm,N,3);
278cac4c232SBarry Smith   PetscTryMethod(dm,"DMRedundantSetSize_C",(DM,PetscMPIInt,PetscInt),(dm,rank,N));
2798ac4e037SJed Brown   PetscFunctionReturn(0);
2808ac4e037SJed Brown }
2818ac4e037SJed Brown 
2828ac4e037SJed Brown /*@
2838ac4e037SJed Brown     DMRedundantGetSize - Gets the size of a densely coupled redundant object
2848ac4e037SJed Brown 
2858ac4e037SJed Brown     Not Collective
2868ac4e037SJed Brown 
2878ac4e037SJed Brown     Input Parameter:
288a2b725a8SWilliam Gropp .   dm - redundant DM
2898ac4e037SJed Brown 
2908ac4e037SJed Brown     Output Parameters:
2910298fd71SBarry Smith +   rank - rank of process to own redundant degrees of freedom (or NULL)
2920298fd71SBarry Smith -   N - total number of redundant degrees of freedom (or NULL)
2938ac4e037SJed Brown 
2948ac4e037SJed Brown     Level: advanced
2958ac4e037SJed Brown 
2968ac4e037SJed Brown .seealso DMDestroy(), DMCreateGlobalVector(), DMRedundantCreate(), DMRedundantSetSize()
2978ac4e037SJed Brown @*/
298907376e6SBarry Smith PetscErrorCode DMRedundantGetSize(DM dm,PetscMPIInt *rank,PetscInt *N)
2998ac4e037SJed Brown {
3008ac4e037SJed Brown   PetscFunctionBegin;
3018ac4e037SJed Brown   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
3028ac4e037SJed Brown   PetscValidType(dm,1);
303cac4c232SBarry Smith   PetscUseMethod(dm,"DMRedundantGetSize_C",(DM,PetscMPIInt*,PetscInt*),(dm,rank,N));
3048ac4e037SJed Brown   PetscFunctionReturn(0);
3058ac4e037SJed Brown }
3068ac4e037SJed Brown 
307907376e6SBarry Smith static PetscErrorCode DMRedundantSetSize_Redundant(DM dm,PetscMPIInt rank,PetscInt N)
3088ac4e037SJed Brown {
3098ac4e037SJed Brown   DM_Redundant   *red = (DM_Redundant*)dm->data;
3108ac4e037SJed Brown   PetscMPIInt    myrank;
31173c9e547SStefano Zampini   PetscInt       i,*globals;
3128ac4e037SJed Brown 
3138ac4e037SJed Brown   PetscFunctionBegin;
3149566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)dm),&myrank));
3158ac4e037SJed Brown   red->rank = rank;
3168ac4e037SJed Brown   red->N    = N;
3178ac4e037SJed Brown   red->n    = (myrank == rank) ? N : 0;
31873c9e547SStefano Zampini 
31973c9e547SStefano Zampini   /* mapping is setup here */
3209566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(red->N,&globals));
32173c9e547SStefano Zampini   for (i=0; i<red->N; i++) globals[i] = i;
3229566063dSJacob Faibussowitsch   PetscCall(ISLocalToGlobalMappingDestroy(&dm->ltogmap));
3239566063dSJacob Faibussowitsch   PetscCall(ISLocalToGlobalMappingCreate(PetscObjectComm((PetscObject)dm),1,red->N,globals,PETSC_OWN_POINTER,&dm->ltogmap));
3248ac4e037SJed Brown   PetscFunctionReturn(0);
3258ac4e037SJed Brown }
3268ac4e037SJed Brown 
3271e6b0712SBarry Smith static PetscErrorCode DMRedundantGetSize_Redundant(DM dm,PetscInt *rank,PetscInt *N)
3288ac4e037SJed Brown {
3298ac4e037SJed Brown   DM_Redundant *red = (DM_Redundant*)dm->data;
3308ac4e037SJed Brown 
3318ac4e037SJed Brown   PetscFunctionBegin;
3328ac4e037SJed Brown   if (rank) *rank = red->rank;
3338ac4e037SJed Brown   if (N)    *N = red->N;
3348ac4e037SJed Brown   PetscFunctionReturn(0);
3358ac4e037SJed Brown }
3368ac4e037SJed Brown 
33736623e74SStefano Zampini static PetscErrorCode DMSetUpGLVisViewer_Redundant(PetscObject odm, PetscViewer viewer)
33836623e74SStefano Zampini {
33936623e74SStefano Zampini   PetscFunctionBegin;
34036623e74SStefano Zampini   PetscFunctionReturn(0);
34136623e74SStefano Zampini }
34236623e74SStefano Zampini 
3433efe6655SBarry Smith /*MC
3443efe6655SBarry Smith    DMREDUNDANT = "redundant" - A DM object that is used to manage data for a small set of dense globally coupled variables.
3453efe6655SBarry Smith          In the global representation of the vector the variables are all stored on a single MPI process (all the other MPI processes have
3463efe6655SBarry Smith          no variables) in the local representation all the variables are stored on ALL the MPI processes (because they are all needed for each
3473efe6655SBarry Smith          processes local computations).
3483efe6655SBarry Smith 
3493efe6655SBarry Smith          This DM is generally used inside a DMCOMPOSITE object. For example, it may be used to store continuation parameters for a bifurcation problem.
3503efe6655SBarry Smith 
3513efe6655SBarry Smith   Level: intermediate
3523efe6655SBarry Smith 
3531abcec8cSBarry Smith .seealso: DMType, DMCOMPOSITE,  DMCreate(), DMRedundantSetSize(), DMRedundantGetSize()
3543efe6655SBarry Smith M*/
3553efe6655SBarry Smith 
3568cc058d9SJed Brown PETSC_EXTERN PetscErrorCode DMCreate_Redundant(DM dm)
3578ac4e037SJed Brown {
3588ac4e037SJed Brown   DM_Redundant   *red;
3598ac4e037SJed Brown 
3608ac4e037SJed Brown   PetscFunctionBegin;
3619566063dSJacob Faibussowitsch   PetscCall(PetscNewLog(dm,&red));
3628ac4e037SJed Brown   dm->data = red;
3638ac4e037SJed Brown 
3648ac4e037SJed Brown   dm->ops->setup               = DMSetUp_Redundant;
3658ac4e037SJed Brown   dm->ops->view                = DMView_Redundant;
3668ac4e037SJed Brown   dm->ops->createglobalvector  = DMCreateGlobalVector_Redundant;
3678ac4e037SJed Brown   dm->ops->createlocalvector   = DMCreateLocalVector_Redundant;
36825296bd5SBarry Smith   dm->ops->creatematrix        = DMCreateMatrix_Redundant;
3698ac4e037SJed Brown   dm->ops->destroy             = DMDestroy_Redundant;
3708ac4e037SJed Brown   dm->ops->globaltolocalbegin  = DMGlobalToLocalBegin_Redundant;
3718ac4e037SJed Brown   dm->ops->globaltolocalend    = DMGlobalToLocalEnd_Redundant;
3728ac4e037SJed Brown   dm->ops->localtoglobalbegin  = DMLocalToGlobalBegin_Redundant;
3738ac4e037SJed Brown   dm->ops->localtoglobalend    = DMLocalToGlobalEnd_Redundant;
3748ac4e037SJed Brown   dm->ops->refine              = DMRefine_Redundant;
3758ac4e037SJed Brown   dm->ops->coarsen             = DMCoarsen_Redundant;
37625296bd5SBarry Smith   dm->ops->createinterpolation = DMCreateInterpolation_Redundant;
377e727c939SJed Brown   dm->ops->getcoloring         = DMCreateColoring_Redundant;
3788865f1eaSKarl Rupp 
3799566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)dm,"DMRedundantSetSize_C",DMRedundantSetSize_Redundant));
3809566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)dm,"DMRedundantGetSize_C",DMRedundantGetSize_Redundant));
3819566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)dm,"DMSetUpGLVisViewer_C",DMSetUpGLVisViewer_Redundant));
3828ac4e037SJed Brown   PetscFunctionReturn(0);
3838ac4e037SJed Brown }
3848ac4e037SJed Brown 
3858ac4e037SJed Brown /*@C
3868ac4e037SJed Brown     DMRedundantCreate - Creates a DM object, used to manage data for dense globally coupled variables
3878ac4e037SJed Brown 
388d083f849SBarry Smith     Collective
3898ac4e037SJed Brown 
390d8d19677SJose E. Roman     Input Parameters:
3918ac4e037SJed Brown +   comm - the processors that will share the global vector
3928ac4e037SJed Brown .   rank - rank to own the redundant values
3938ac4e037SJed Brown -   N - total number of degrees of freedom
3948ac4e037SJed Brown 
3958ac4e037SJed Brown     Output Parameters:
396907376e6SBarry Smith .   dm - the redundant DM
3978ac4e037SJed Brown 
3988ac4e037SJed Brown     Level: advanced
3998ac4e037SJed Brown 
4003efe6655SBarry Smith .seealso DMDestroy(), DMCreateGlobalVector(), DMCreateMatrix(), DMCompositeAddDM(), DMREDUNDANT, DMSetType(), DMRedundantSetSize(), DMRedundantGetSize()
4018ac4e037SJed Brown 
4028ac4e037SJed Brown @*/
403907376e6SBarry Smith PetscErrorCode DMRedundantCreate(MPI_Comm comm,PetscMPIInt rank,PetscInt N,DM *dm)
4048ac4e037SJed Brown {
4058ac4e037SJed Brown   PetscFunctionBegin;
406064a246eSJacob Faibussowitsch   PetscValidPointer(dm,4);
4079566063dSJacob Faibussowitsch   PetscCall(DMCreate(comm,dm));
4089566063dSJacob Faibussowitsch   PetscCall(DMSetType(*dm,DMREDUNDANT));
4099566063dSJacob Faibussowitsch   PetscCall(DMRedundantSetSize(*dm,rank,N));
4109566063dSJacob Faibussowitsch   PetscCall(DMSetUp(*dm));
4118ac4e037SJed Brown   PetscFunctionReturn(0);
4128ac4e037SJed Brown }
413