xref: /petsc/src/dm/impls/redundant/dmredundant.c (revision 45b6f7e9ecd564e8bb535880b270d4a8084ff211)
107475bc1SBarry Smith #include <petsc-private/dmimpl.h>
28ac4e037SJed Brown #include <petscdmredundant.h>   /*I      "petscdmredundant.h" I*/
38ac4e037SJed Brown 
48ac4e037SJed Brown typedef struct  {
58ac4e037SJed Brown   PetscInt 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 
108ac4e037SJed Brown #undef __FUNCT__
11950540a4SJed Brown #define __FUNCT__ "DMCreateMatrix_Redundant"
12b412c318SBarry Smith static PetscErrorCode DMCreateMatrix_Redundant(DM dm,Mat *J)
138ac4e037SJed Brown {
148ac4e037SJed Brown   DM_Redundant           *red = (DM_Redundant*)dm->data;
158ac4e037SJed Brown   PetscErrorCode         ierr;
16*45b6f7e9SBarry Smith   ISLocalToGlobalMapping ltog;
17d2e2b695SJed Brown   PetscInt               i,rstart,rend,*cols;
18d2e2b695SJed Brown   PetscScalar            *vals;
198ac4e037SJed Brown 
208ac4e037SJed Brown   PetscFunctionBegin;
21ce94432eSBarry Smith   ierr = MatCreate(PetscObjectComm((PetscObject)dm),J);CHKERRQ(ierr);
228ac4e037SJed Brown   ierr = MatSetSizes(*J,red->n,red->n,red->N,red->N);CHKERRQ(ierr);
23b412c318SBarry Smith   ierr = MatSetType(*J,dm->mattype);CHKERRQ(ierr);
240298fd71SBarry Smith   ierr = MatSeqAIJSetPreallocation(*J,red->n,NULL);CHKERRQ(ierr);
250298fd71SBarry Smith   ierr = MatSeqBAIJSetPreallocation(*J,1,red->n,NULL);CHKERRQ(ierr);
260298fd71SBarry Smith   ierr = MatMPIAIJSetPreallocation(*J,red->n,NULL,red->N-red->n,NULL);CHKERRQ(ierr);
270298fd71SBarry Smith   ierr = MatMPIBAIJSetPreallocation(*J,1,red->n,NULL,red->N-red->n,NULL);CHKERRQ(ierr);
288ac4e037SJed Brown 
298ac4e037SJed Brown   ierr = DMGetLocalToGlobalMapping(dm,&ltog);CHKERRQ(ierr);
308ac4e037SJed Brown   ierr = MatSetLocalToGlobalMapping(*J,ltog,ltog);CHKERRQ(ierr);
31d2e2b695SJed Brown 
32dcca6d9dSJed Brown   ierr = PetscMalloc2(red->N,&cols,red->N,&vals);CHKERRQ(ierr);
33d2e2b695SJed Brown   for (i=0; i<red->N; i++) {
34d2e2b695SJed Brown     cols[i] = i;
35d2e2b695SJed Brown     vals[i] = 0.0;
36d2e2b695SJed Brown   }
37d2e2b695SJed Brown   ierr = MatGetOwnershipRange(*J,&rstart,&rend);CHKERRQ(ierr);
38d2e2b695SJed Brown   for (i=rstart; i<rend; i++) {
39d2e2b695SJed Brown     ierr = MatSetValues(*J,1,&i,red->N,cols,vals,INSERT_VALUES);CHKERRQ(ierr);
40d2e2b695SJed Brown   }
41d2e2b695SJed Brown   ierr = PetscFree2(cols,vals);CHKERRQ(ierr);
42d2e2b695SJed Brown   ierr = MatAssemblyBegin(*J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
43d2e2b695SJed Brown   ierr = MatAssemblyEnd(*J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
448ac4e037SJed Brown   PetscFunctionReturn(0);
458ac4e037SJed Brown }
468ac4e037SJed Brown 
478ac4e037SJed Brown #undef __FUNCT__
488ac4e037SJed Brown #define __FUNCT__ "DMDestroy_Redundant"
498ac4e037SJed Brown static PetscErrorCode DMDestroy_Redundant(DM dm)
508ac4e037SJed Brown {
518ac4e037SJed Brown   PetscErrorCode ierr;
528ac4e037SJed Brown 
538ac4e037SJed Brown   PetscFunctionBegin;
54bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)dm,"DMRedundantSetSize_C",NULL);CHKERRQ(ierr);
55bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)dm,"DMRedundantGetSize_C",NULL);CHKERRQ(ierr);
56435a35e8SMatthew G Knepley   /* This was originally freed in DMDestroy(), but that prevents reference counting of backend objects */
57435a35e8SMatthew G Knepley   ierr = PetscFree(dm->data);CHKERRQ(ierr);
588ac4e037SJed Brown   PetscFunctionReturn(0);
598ac4e037SJed Brown }
608ac4e037SJed Brown 
618ac4e037SJed Brown #undef __FUNCT__
628ac4e037SJed Brown #define __FUNCT__ "DMCreateGlobalVector_Redundant"
638ac4e037SJed Brown static PetscErrorCode DMCreateGlobalVector_Redundant(DM dm,Vec *gvec)
648ac4e037SJed Brown {
658ac4e037SJed Brown   PetscErrorCode         ierr;
668ac4e037SJed Brown   DM_Redundant           *red = (DM_Redundant*)dm->data;
67*45b6f7e9SBarry Smith   ISLocalToGlobalMapping ltog;
688ac4e037SJed Brown 
698ac4e037SJed Brown   PetscFunctionBegin;
708ac4e037SJed Brown   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
718ac4e037SJed Brown   PetscValidPointer(gvec,2);
728ac4e037SJed Brown   *gvec = 0;
73ce94432eSBarry Smith   ierr  = VecCreate(PetscObjectComm((PetscObject)dm),gvec);CHKERRQ(ierr);
748ac4e037SJed Brown   ierr  = VecSetSizes(*gvec,red->n,red->N);CHKERRQ(ierr);
758ac4e037SJed Brown   ierr  = VecSetType(*gvec,dm->vectype);CHKERRQ(ierr);
768ac4e037SJed Brown   ierr  = DMGetLocalToGlobalMapping(dm,&ltog);CHKERRQ(ierr);
778ac4e037SJed Brown   ierr  = VecSetLocalToGlobalMapping(*gvec,ltog);CHKERRQ(ierr);
78c688c046SMatthew G Knepley   ierr  = VecSetDM(*gvec,dm);CHKERRQ(ierr);
798ac4e037SJed Brown   PetscFunctionReturn(0);
808ac4e037SJed Brown }
818ac4e037SJed Brown 
828ac4e037SJed Brown #undef __FUNCT__
838ac4e037SJed Brown #define __FUNCT__ "DMCreateLocalVector_Redundant"
848ac4e037SJed Brown static PetscErrorCode DMCreateLocalVector_Redundant(DM dm,Vec *lvec)
858ac4e037SJed Brown {
868ac4e037SJed Brown   PetscErrorCode ierr;
878ac4e037SJed Brown   DM_Redundant   *red = (DM_Redundant*)dm->data;
888ac4e037SJed Brown 
898ac4e037SJed Brown   PetscFunctionBegin;
908ac4e037SJed Brown   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
918ac4e037SJed Brown   PetscValidPointer(lvec,2);
928ac4e037SJed Brown   *lvec = 0;
938ac4e037SJed Brown   ierr  = VecCreate(PETSC_COMM_SELF,lvec);CHKERRQ(ierr);
948ac4e037SJed Brown   ierr  = VecSetSizes(*lvec,red->N,red->N);CHKERRQ(ierr);
958ac4e037SJed Brown   ierr  = VecSetType(*lvec,dm->vectype);CHKERRQ(ierr);
96c688c046SMatthew G Knepley   ierr  = VecSetDM(*lvec,dm);CHKERRQ(ierr);
978ac4e037SJed Brown   PetscFunctionReturn(0);
988ac4e037SJed Brown }
998ac4e037SJed Brown 
1008ac4e037SJed Brown #undef __FUNCT__
1018ac4e037SJed Brown #define __FUNCT__ "DMLocalToGlobalBegin_Redundant"
1028ac4e037SJed Brown static PetscErrorCode DMLocalToGlobalBegin_Redundant(DM dm,Vec l,InsertMode imode,Vec g)
1038ac4e037SJed Brown {
1048ac4e037SJed Brown   PetscErrorCode    ierr;
1058ac4e037SJed Brown   DM_Redundant      *red = (DM_Redundant*)dm->data;
1068ac4e037SJed Brown   const PetscScalar *lv;
1078ac4e037SJed Brown   PetscScalar       *gv;
1088ac4e037SJed Brown   PetscMPIInt       rank;
1098ac4e037SJed Brown 
1108ac4e037SJed Brown   PetscFunctionBegin;
111ce94432eSBarry Smith   ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)dm),&rank);CHKERRQ(ierr);
1128ac4e037SJed Brown   ierr = VecGetArrayRead(l,&lv);CHKERRQ(ierr);
1138ac4e037SJed Brown   ierr = VecGetArray(g,&gv);CHKERRQ(ierr);
1148ac4e037SJed Brown   switch (imode) {
1158ac4e037SJed Brown   case ADD_VALUES:
1168ac4e037SJed Brown   case MAX_VALUES:
1178ac4e037SJed Brown   {
1188ac4e037SJed Brown     void        *source;
1198ac4e037SJed Brown     PetscScalar *buffer;
1208ac4e037SJed Brown     PetscInt    i;
1218ac4e037SJed Brown     if (rank == red->rank) {
1228ac4e037SJed Brown #if defined(PETSC_HAVE_MPI_IN_PLACE)
1238ac4e037SJed Brown       buffer = gv;
1248ac4e037SJed Brown       source = MPI_IN_PLACE;
1258ac4e037SJed Brown #else
126785e854fSJed Brown       ierr   = PetscMalloc1(red->N,&buffer);CHKERRQ(ierr);
1278ac4e037SJed Brown       source = buffer;
1288ac4e037SJed Brown #endif
1298ac4e037SJed Brown       if (imode == ADD_VALUES) for (i=0; i<red->N; i++) buffer[i] = gv[i] + lv[i];
1308ac4e037SJed Brown #if !defined(PETSC_USE_COMPLEX)
1318ac4e037SJed Brown       if (imode == MAX_VALUES) for (i=0; i<red->N; i++) buffer[i] = PetscMax(gv[i],lv[i]);
1328ac4e037SJed Brown #endif
1338865f1eaSKarl Rupp     } else source = (void*)lv;
1348fa295daSBarry Smith     ierr = MPI_Reduce(source,gv,red->N,MPIU_SCALAR,(imode == ADD_VALUES) ? MPIU_SUM : MPIU_MAX,red->rank,PetscObjectComm((PetscObject)dm));CHKERRQ(ierr);
1358ac4e037SJed Brown #if !defined(PETSC_HAVE_MPI_IN_PLACE)
13699839c1eSBarry Smith     if (rank == red->rank) {ierr = PetscFree(buffer);CHKERRQ(ierr);}
1378ac4e037SJed Brown #endif
1388ac4e037SJed Brown   } break;
1398ac4e037SJed Brown   case INSERT_VALUES:
1408ac4e037SJed Brown     ierr = PetscMemcpy(gv,lv,red->n*sizeof(PetscScalar));CHKERRQ(ierr);
1418ac4e037SJed Brown     break;
142ce94432eSBarry Smith   default: SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"InsertMode not supported");
1438ac4e037SJed Brown   }
1448ac4e037SJed Brown   ierr = VecRestoreArrayRead(l,&lv);CHKERRQ(ierr);
1458ac4e037SJed Brown   ierr = VecRestoreArray(g,&gv);CHKERRQ(ierr);
1468ac4e037SJed Brown   PetscFunctionReturn(0);
1478ac4e037SJed Brown }
1488ac4e037SJed Brown 
1498ac4e037SJed Brown #undef __FUNCT__
1508ac4e037SJed Brown #define __FUNCT__ "DMLocalToGlobalEnd_Redundant"
1518ac4e037SJed Brown static PetscErrorCode DMLocalToGlobalEnd_Redundant(DM dm,Vec l,InsertMode imode,Vec g)
1528ac4e037SJed Brown {
1538ac4e037SJed Brown   PetscFunctionBegin;
1548ac4e037SJed Brown   PetscFunctionReturn(0);
1558ac4e037SJed Brown }
1568ac4e037SJed Brown 
1578ac4e037SJed Brown #undef __FUNCT__
1588ac4e037SJed Brown #define __FUNCT__ "DMGlobalToLocalBegin_Redundant"
1598ac4e037SJed Brown static PetscErrorCode DMGlobalToLocalBegin_Redundant(DM dm,Vec g,InsertMode imode,Vec l)
1608ac4e037SJed Brown {
1618ac4e037SJed Brown   PetscErrorCode    ierr;
1628ac4e037SJed Brown   DM_Redundant      *red = (DM_Redundant*)dm->data;
1638ac4e037SJed Brown   const PetscScalar *gv;
1648ac4e037SJed Brown   PetscScalar       *lv;
1658ac4e037SJed Brown 
1668ac4e037SJed Brown   PetscFunctionBegin;
1678ac4e037SJed Brown   ierr = VecGetArrayRead(g,&gv);CHKERRQ(ierr);
1688ac4e037SJed Brown   ierr = VecGetArray(l,&lv);CHKERRQ(ierr);
1698ac4e037SJed Brown   switch (imode) {
1708ac4e037SJed Brown   case INSERT_VALUES:
1718ac4e037SJed Brown     if (red->n) {ierr = PetscMemcpy(lv,gv,red->n*sizeof(PetscScalar));CHKERRQ(ierr);}
172ce94432eSBarry Smith     ierr = MPI_Bcast(lv,red->N,MPIU_SCALAR,red->rank,PetscObjectComm((PetscObject)dm));CHKERRQ(ierr);
1738ac4e037SJed Brown     break;
174ce94432eSBarry Smith   default: SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"InsertMode not supported");
1758ac4e037SJed Brown   }
1768ac4e037SJed Brown   ierr = VecRestoreArrayRead(g,&gv);CHKERRQ(ierr);
1778ac4e037SJed Brown   ierr = VecRestoreArray(l,&lv);CHKERRQ(ierr);
1788ac4e037SJed Brown   PetscFunctionReturn(0);
1798ac4e037SJed Brown }
1808ac4e037SJed Brown 
1818ac4e037SJed Brown #undef __FUNCT__
1828ac4e037SJed Brown #define __FUNCT__ "DMGlobalToLocalEnd_Redundant"
1838ac4e037SJed Brown static PetscErrorCode DMGlobalToLocalEnd_Redundant(DM dm,Vec g,InsertMode imode,Vec l)
1848ac4e037SJed Brown {
1858ac4e037SJed Brown   PetscFunctionBegin;
1868ac4e037SJed Brown   PetscFunctionReturn(0);
1878ac4e037SJed Brown }
1888ac4e037SJed Brown 
1898ac4e037SJed Brown #undef __FUNCT__
1908ac4e037SJed Brown #define __FUNCT__ "DMSetUp_Redundant"
1918ac4e037SJed Brown static PetscErrorCode DMSetUp_Redundant(DM dm)
1928ac4e037SJed Brown {
1938ac4e037SJed Brown   PetscErrorCode ierr;
1948ac4e037SJed Brown   DM_Redundant   *red = (DM_Redundant*)dm->data;
1958ac4e037SJed Brown   PetscInt       i,*globals;
1968ac4e037SJed Brown 
1978ac4e037SJed Brown   PetscFunctionBegin;
198785e854fSJed Brown   ierr = PetscMalloc1(red->N,&globals);CHKERRQ(ierr);
1998ac4e037SJed Brown   for (i=0; i<red->N; i++) globals[i] = i;
200f0413b6fSBarry Smith   ierr         = ISLocalToGlobalMappingCreate(PETSC_COMM_SELF,1,red->N,globals,PETSC_OWN_POINTER,&dm->ltogmap);CHKERRQ(ierr);
2018ac4e037SJed Brown   ierr         = PetscObjectReference((PetscObject)dm->ltogmap);CHKERRQ(ierr);
2028ac4e037SJed Brown   PetscFunctionReturn(0);
2038ac4e037SJed Brown }
2048ac4e037SJed Brown 
2058ac4e037SJed Brown #undef __FUNCT__
2068ac4e037SJed Brown #define __FUNCT__ "DMView_Redundant"
2078ac4e037SJed Brown static PetscErrorCode DMView_Redundant(DM dm,PetscViewer viewer)
2088ac4e037SJed Brown {
2098ac4e037SJed Brown   PetscErrorCode ierr;
2108ac4e037SJed Brown   DM_Redundant   *red = (DM_Redundant*)dm->data;
2118ac4e037SJed Brown   PetscBool      iascii;
2128ac4e037SJed Brown 
2138ac4e037SJed Brown   PetscFunctionBegin;
214251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
2158ac4e037SJed Brown   if (iascii) {
2168ac4e037SJed Brown     ierr = PetscViewerASCIIPrintf(viewer,"redundant: rank=%D N=%D\n",red->rank,red->N);CHKERRQ(ierr);
2178ac4e037SJed Brown   }
2188ac4e037SJed Brown   PetscFunctionReturn(0);
2198ac4e037SJed Brown }
2208ac4e037SJed Brown 
22169527181SJed Brown #undef __FUNCT__
222e727c939SJed Brown #define __FUNCT__ "DMCreateColoring_Redundant"
223b412c318SBarry Smith static PetscErrorCode DMCreateColoring_Redundant(DM dm,ISColoringType ctype,ISColoring *coloring)
22469527181SJed Brown {
22569527181SJed Brown   DM_Redundant    *red = (DM_Redundant*)dm->data;
22669527181SJed Brown   PetscErrorCode  ierr;
22769527181SJed Brown   PetscInt        i,nloc;
22869527181SJed Brown   ISColoringValue *colors;
22969527181SJed Brown 
23069527181SJed Brown   PetscFunctionBegin;
23169527181SJed Brown   switch (ctype) {
23269527181SJed Brown   case IS_COLORING_GLOBAL:
23369527181SJed Brown     nloc = red->n;
23469527181SJed Brown     break;
23569527181SJed Brown   case IS_COLORING_GHOSTED:
23669527181SJed Brown     nloc = red->N;
23769527181SJed Brown     break;
238ce94432eSBarry Smith   default: SETERRQ1(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONG,"Unknown ISColoringType %d",(int)ctype);
23969527181SJed Brown   }
240785e854fSJed Brown   ierr = PetscMalloc1(nloc,&colors);CHKERRQ(ierr);
24169527181SJed Brown   for (i=0; i<nloc; i++) colors[i] = i;
242ce94432eSBarry Smith   ierr = ISColoringCreate(PetscObjectComm((PetscObject)dm),red->N,nloc,colors,coloring);CHKERRQ(ierr);
24369527181SJed Brown   ierr = ISColoringSetType(*coloring,ctype);CHKERRQ(ierr);
24469527181SJed Brown   PetscFunctionReturn(0);
24569527181SJed Brown }
2468ac4e037SJed Brown 
2478ac4e037SJed Brown #undef __FUNCT__
2488ac4e037SJed Brown #define __FUNCT__ "DMRefine_Redundant"
2498ac4e037SJed Brown static PetscErrorCode DMRefine_Redundant(DM dmc,MPI_Comm comm,DM *dmf)
2508ac4e037SJed Brown {
2518ac4e037SJed Brown   PetscErrorCode ierr;
2528ac4e037SJed Brown   PetscMPIInt    flag;
2538ac4e037SJed Brown   DM_Redundant   *redc = (DM_Redundant*)dmc->data;
2548ac4e037SJed Brown 
2558ac4e037SJed Brown   PetscFunctionBegin;
256ce94432eSBarry Smith   if (comm == MPI_COMM_NULL) {
257ce94432eSBarry Smith     ierr = PetscObjectGetComm((PetscObject)dmc,&comm);CHKERRQ(ierr);
258ce94432eSBarry Smith   }
259ce94432eSBarry Smith   ierr = MPI_Comm_compare(PetscObjectComm((PetscObject)dmc),comm,&flag);CHKERRQ(ierr);
260ce94432eSBarry Smith   if (flag != MPI_CONGRUENT && flag != MPI_IDENT) SETERRQ(PetscObjectComm((PetscObject)dmc),PETSC_ERR_SUP,"cannot change communicators");
2618ac4e037SJed Brown   ierr = DMRedundantCreate(comm,redc->rank,redc->N,dmf);CHKERRQ(ierr);
2628ac4e037SJed Brown   PetscFunctionReturn(0);
2638ac4e037SJed Brown }
2648ac4e037SJed Brown 
2658ac4e037SJed Brown #undef __FUNCT__
2668ac4e037SJed Brown #define __FUNCT__ "DMCoarsen_Redundant"
2678ac4e037SJed Brown static PetscErrorCode DMCoarsen_Redundant(DM dmf,MPI_Comm comm,DM *dmc)
2688ac4e037SJed Brown {
2698ac4e037SJed Brown   PetscErrorCode ierr;
2708ac4e037SJed Brown   PetscMPIInt    flag;
2718ac4e037SJed Brown   DM_Redundant   *redf = (DM_Redundant*)dmf->data;
2728ac4e037SJed Brown 
2738ac4e037SJed Brown   PetscFunctionBegin;
274ce94432eSBarry Smith   if (comm == MPI_COMM_NULL) {
275ce94432eSBarry Smith     ierr = PetscObjectGetComm((PetscObject)dmf,&comm);CHKERRQ(ierr);
276ce94432eSBarry Smith   }
277ce94432eSBarry Smith   ierr = MPI_Comm_compare(PetscObjectComm((PetscObject)dmf),comm,&flag);CHKERRQ(ierr);
278ce94432eSBarry Smith   if (flag != MPI_CONGRUENT && flag != MPI_IDENT) SETERRQ(PetscObjectComm((PetscObject)dmf),PETSC_ERR_SUP,"cannot change communicators");
2798ac4e037SJed Brown   ierr = DMRedundantCreate(comm,redf->rank,redf->N,dmc);CHKERRQ(ierr);
2808ac4e037SJed Brown   PetscFunctionReturn(0);
2818ac4e037SJed Brown }
2828ac4e037SJed Brown 
2838ac4e037SJed Brown #undef __FUNCT__
284e727c939SJed Brown #define __FUNCT__ "DMCreateInterpolation_Redundant"
285e727c939SJed Brown static PetscErrorCode DMCreateInterpolation_Redundant(DM dmc,DM dmf,Mat *P,Vec *scale)
2868ac4e037SJed Brown {
2878ac4e037SJed Brown   PetscErrorCode ierr;
2888ac4e037SJed Brown   DM_Redundant   *redc = (DM_Redundant*)dmc->data;
2898ac4e037SJed Brown   DM_Redundant   *redf = (DM_Redundant*)dmf->data;
2908ac4e037SJed Brown   PetscMPIInt    flag;
2918ac4e037SJed Brown   PetscInt       i,rstart,rend;
2928ac4e037SJed Brown 
2938ac4e037SJed Brown   PetscFunctionBegin;
294ce94432eSBarry Smith   ierr = MPI_Comm_compare(PetscObjectComm((PetscObject)dmc),PetscObjectComm((PetscObject)dmf),&flag);CHKERRQ(ierr);
295ce94432eSBarry Smith   if (flag != MPI_CONGRUENT && flag != MPI_IDENT) SETERRQ(PetscObjectComm((PetscObject)dmf),PETSC_ERR_SUP,"cannot change communicators");
296ce94432eSBarry Smith   if (redc->rank != redf->rank) SETERRQ(PetscObjectComm((PetscObject)dmf),PETSC_ERR_ARG_INCOMP,"Owning rank does not match");
297ce94432eSBarry Smith   if (redc->N != redf->N) SETERRQ(PetscObjectComm((PetscObject)dmf),PETSC_ERR_ARG_INCOMP,"Global size does not match");
298ce94432eSBarry Smith   ierr = MatCreate(PetscObjectComm((PetscObject)dmc),P);CHKERRQ(ierr);
2998ac4e037SJed Brown   ierr = MatSetSizes(*P,redc->n,redc->n,redc->N,redc->N);CHKERRQ(ierr);
3008ac4e037SJed Brown   ierr = MatSetType(*P,MATAIJ);CHKERRQ(ierr);
3018ac4e037SJed Brown   ierr = MatSeqAIJSetPreallocation(*P,1,0);CHKERRQ(ierr);
3028ac4e037SJed Brown   ierr = MatMPIAIJSetPreallocation(*P,1,0,0,0);CHKERRQ(ierr);
3038ac4e037SJed Brown   ierr = MatGetOwnershipRange(*P,&rstart,&rend);CHKERRQ(ierr);
3048ac4e037SJed Brown   for (i=rstart; i<rend; i++) {ierr = MatSetValue(*P,i,i,1.0,INSERT_VALUES);CHKERRQ(ierr);}
3058ac4e037SJed Brown   ierr = MatAssemblyBegin(*P,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
3068ac4e037SJed Brown   ierr = MatAssemblyEnd(*P,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
307e727c939SJed Brown   if (scale) {ierr = DMCreateInterpolationScale(dmc,dmf,*P,scale);CHKERRQ(ierr);}
3088ac4e037SJed Brown   PetscFunctionReturn(0);
3098ac4e037SJed Brown }
3108ac4e037SJed Brown 
3118ac4e037SJed Brown #undef __FUNCT__
3128ac4e037SJed Brown #define __FUNCT__ "DMRedundantSetSize"
3138ac4e037SJed Brown /*@
3148ac4e037SJed Brown     DMRedundantSetSize - Sets the size of a densely coupled redundant object
3158ac4e037SJed Brown 
3168ac4e037SJed Brown     Collective on DM
3178ac4e037SJed Brown 
3188ac4e037SJed Brown     Input Parameter:
3198ac4e037SJed Brown +   dm - redundant DM
3208ac4e037SJed Brown .   rank - rank of process to own redundant degrees of freedom
3218ac4e037SJed Brown -   N - total number of redundant degrees of freedom
3228ac4e037SJed Brown 
3238ac4e037SJed Brown     Level: advanced
3248ac4e037SJed Brown 
3258ac4e037SJed Brown .seealso DMDestroy(), DMCreateGlobalVector(), DMRedundantCreate(), DMRedundantGetSize()
3268ac4e037SJed Brown @*/
3278ac4e037SJed Brown PetscErrorCode DMRedundantSetSize(DM dm,PetscInt rank,PetscInt N)
3288ac4e037SJed Brown {
3298ac4e037SJed Brown   PetscErrorCode ierr;
3308ac4e037SJed Brown 
3318ac4e037SJed Brown   PetscFunctionBegin;
3328ac4e037SJed Brown   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
3338ac4e037SJed Brown   PetscValidType(dm,1);
3348ac4e037SJed Brown   PetscValidLogicalCollectiveInt(dm,rank,2);
3358ac4e037SJed Brown   PetscValidLogicalCollectiveInt(dm,N,3);
3368ac4e037SJed Brown   ierr = PetscTryMethod(dm,"DMRedundantSetSize_C",(DM,PetscInt,PetscInt),(dm,rank,N));CHKERRQ(ierr);
3378ac4e037SJed Brown   PetscFunctionReturn(0);
3388ac4e037SJed Brown }
3398ac4e037SJed Brown 
3408ac4e037SJed Brown #undef __FUNCT__
3418ac4e037SJed Brown #define __FUNCT__ "DMRedundantGetSize"
3428ac4e037SJed Brown /*@
3438ac4e037SJed Brown     DMRedundantGetSize - Gets the size of a densely coupled redundant object
3448ac4e037SJed Brown 
3458ac4e037SJed Brown     Not Collective
3468ac4e037SJed Brown 
3478ac4e037SJed Brown     Input Parameter:
3488ac4e037SJed Brown +   dm - redundant DM
3498ac4e037SJed Brown 
3508ac4e037SJed Brown     Output Parameters:
3510298fd71SBarry Smith +   rank - rank of process to own redundant degrees of freedom (or NULL)
3520298fd71SBarry Smith -   N - total number of redundant degrees of freedom (or NULL)
3538ac4e037SJed Brown 
3548ac4e037SJed Brown     Level: advanced
3558ac4e037SJed Brown 
3568ac4e037SJed Brown .seealso DMDestroy(), DMCreateGlobalVector(), DMRedundantCreate(), DMRedundantSetSize()
3578ac4e037SJed Brown @*/
3588ac4e037SJed Brown PetscErrorCode DMRedundantGetSize(DM dm,PetscInt *rank,PetscInt *N)
3598ac4e037SJed Brown {
3608ac4e037SJed Brown   PetscErrorCode ierr;
3618ac4e037SJed Brown 
3628ac4e037SJed Brown   PetscFunctionBegin;
3638ac4e037SJed Brown   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
3648ac4e037SJed Brown   PetscValidType(dm,1);
3658ac4e037SJed Brown   ierr = PetscUseMethod(dm,"DMRedundantGetSize_C",(DM,PetscInt*,PetscInt*),(dm,rank,N));CHKERRQ(ierr);
3668ac4e037SJed Brown   PetscFunctionReturn(0);
3678ac4e037SJed Brown }
3688ac4e037SJed Brown 
3698ac4e037SJed Brown #undef __FUNCT__
3708ac4e037SJed Brown #define __FUNCT__ "DMRedundantSetSize_Redundant"
3711e6b0712SBarry Smith static PetscErrorCode DMRedundantSetSize_Redundant(DM dm,PetscInt rank,PetscInt N)
3728ac4e037SJed Brown {
3738ac4e037SJed Brown   DM_Redundant   *red = (DM_Redundant*)dm->data;
3748ac4e037SJed Brown   PetscErrorCode ierr;
3758ac4e037SJed Brown   PetscMPIInt    myrank;
3768ac4e037SJed Brown 
3778ac4e037SJed Brown   PetscFunctionBegin;
378ce94432eSBarry Smith   ierr      = MPI_Comm_rank(PetscObjectComm((PetscObject)dm),&myrank);CHKERRQ(ierr);
3798ac4e037SJed Brown   red->rank = rank;
3808ac4e037SJed Brown   red->N    = N;
3818ac4e037SJed Brown   red->n    = (myrank == rank) ? N : 0;
3828ac4e037SJed Brown   PetscFunctionReturn(0);
3838ac4e037SJed Brown }
3848ac4e037SJed Brown 
3858ac4e037SJed Brown #undef __FUNCT__
3868ac4e037SJed Brown #define __FUNCT__ "DMRedundantGetSize_Redundant"
3871e6b0712SBarry Smith static PetscErrorCode DMRedundantGetSize_Redundant(DM dm,PetscInt *rank,PetscInt *N)
3888ac4e037SJed Brown {
3898ac4e037SJed Brown   DM_Redundant *red = (DM_Redundant*)dm->data;
3908ac4e037SJed Brown 
3918ac4e037SJed Brown   PetscFunctionBegin;
3928ac4e037SJed Brown   if (rank) *rank = red->rank;
3938ac4e037SJed Brown   if (N)    *N = red->N;
3948ac4e037SJed Brown   PetscFunctionReturn(0);
3958ac4e037SJed Brown }
3968ac4e037SJed Brown 
3973efe6655SBarry Smith /*MC
3983efe6655SBarry Smith    DMREDUNDANT = "redundant" - A DM object that is used to manage data for a small set of dense globally coupled variables.
3993efe6655SBarry Smith          In the global representation of the vector the variables are all stored on a single MPI process (all the other MPI processes have
4003efe6655SBarry Smith          no variables) in the local representation all the variables are stored on ALL the MPI processes (because they are all needed for each
4013efe6655SBarry Smith          processes local computations).
4023efe6655SBarry Smith 
4033efe6655SBarry Smith          This DM is generally used inside a DMCOMPOSITE object. For example, it may be used to store continuation parameters for a bifurcation problem.
4043efe6655SBarry Smith 
4053efe6655SBarry Smith 
4063efe6655SBarry Smith   Level: intermediate
4073efe6655SBarry Smith 
4083efe6655SBarry Smith .seealso: DMType, DMCOMPOSITE, DMCreateRedundant(), DMCreate(), DMRedundantSetSize(), DMRedundantGetSize()
4093efe6655SBarry Smith M*/
4103efe6655SBarry Smith 
4118ac4e037SJed Brown #undef __FUNCT__
4128ac4e037SJed Brown #define __FUNCT__ "DMCreate_Redundant"
4138cc058d9SJed Brown PETSC_EXTERN PetscErrorCode DMCreate_Redundant(DM dm)
4148ac4e037SJed Brown {
4158ac4e037SJed Brown   PetscErrorCode ierr;
4168ac4e037SJed Brown   DM_Redundant   *red;
4178ac4e037SJed Brown 
4188ac4e037SJed Brown   PetscFunctionBegin;
419b00a9115SJed Brown   ierr     = PetscNewLog(dm,&red);CHKERRQ(ierr);
4208ac4e037SJed Brown   dm->data = red;
4218ac4e037SJed Brown 
4228ac4e037SJed Brown   ierr = PetscObjectChangeTypeName((PetscObject)dm,DMREDUNDANT);CHKERRQ(ierr);
4238865f1eaSKarl Rupp 
4248ac4e037SJed Brown   dm->ops->setup              = DMSetUp_Redundant;
4258ac4e037SJed Brown   dm->ops->view               = DMView_Redundant;
4268ac4e037SJed Brown   dm->ops->createglobalvector = DMCreateGlobalVector_Redundant;
4278ac4e037SJed Brown   dm->ops->createlocalvector  = DMCreateLocalVector_Redundant;
42825296bd5SBarry Smith   dm->ops->creatematrix       = DMCreateMatrix_Redundant;
4298ac4e037SJed Brown   dm->ops->destroy            = DMDestroy_Redundant;
4308ac4e037SJed Brown   dm->ops->globaltolocalbegin = DMGlobalToLocalBegin_Redundant;
4318ac4e037SJed Brown   dm->ops->globaltolocalend   = DMGlobalToLocalEnd_Redundant;
4328ac4e037SJed Brown   dm->ops->localtoglobalbegin = DMLocalToGlobalBegin_Redundant;
4338ac4e037SJed Brown   dm->ops->localtoglobalend   = DMLocalToGlobalEnd_Redundant;
4348ac4e037SJed Brown   dm->ops->refine             = DMRefine_Redundant;
4358ac4e037SJed Brown   dm->ops->coarsen            = DMCoarsen_Redundant;
43625296bd5SBarry Smith   dm->ops->createinterpolation= DMCreateInterpolation_Redundant;
437e727c939SJed Brown   dm->ops->getcoloring        = DMCreateColoring_Redundant;
4388865f1eaSKarl Rupp 
439bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)dm,"DMRedundantSetSize_C",DMRedundantSetSize_Redundant);CHKERRQ(ierr);
440bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)dm,"DMRedundantGetSize_C",DMRedundantGetSize_Redundant);CHKERRQ(ierr);
4418ac4e037SJed Brown   PetscFunctionReturn(0);
4428ac4e037SJed Brown }
4438ac4e037SJed Brown 
4448ac4e037SJed Brown #undef __FUNCT__
4458ac4e037SJed Brown #define __FUNCT__ "DMRedundantCreate"
4468ac4e037SJed Brown /*@C
4478ac4e037SJed Brown     DMRedundantCreate - Creates a DM object, used to manage data for dense globally coupled variables
4488ac4e037SJed Brown 
4498ac4e037SJed Brown     Collective on MPI_Comm
4508ac4e037SJed Brown 
4518ac4e037SJed Brown     Input Parameter:
4528ac4e037SJed Brown +   comm - the processors that will share the global vector
4538ac4e037SJed Brown .   rank - rank to own the redundant values
4548ac4e037SJed Brown -   N - total number of degrees of freedom
4558ac4e037SJed Brown 
4568ac4e037SJed Brown     Output Parameters:
4578ac4e037SJed Brown .   red - the redundant DM
4588ac4e037SJed Brown 
4598ac4e037SJed Brown     Level: advanced
4608ac4e037SJed Brown 
4613efe6655SBarry Smith .seealso DMDestroy(), DMCreateGlobalVector(), DMCreateMatrix(), DMCompositeAddDM(), DMREDUNDANT, DMSetType(), DMRedundantSetSize(), DMRedundantGetSize()
4628ac4e037SJed Brown 
4638ac4e037SJed Brown @*/
4648ac4e037SJed Brown PetscErrorCode DMRedundantCreate(MPI_Comm comm,PetscInt rank,PetscInt N,DM *dm)
4658ac4e037SJed Brown {
4668ac4e037SJed Brown   PetscErrorCode ierr;
4678ac4e037SJed Brown 
4688ac4e037SJed Brown   PetscFunctionBegin;
4698ac4e037SJed Brown   PetscValidPointer(dm,2);
4708ac4e037SJed Brown   ierr = DMCreate(comm,dm);CHKERRQ(ierr);
4718ac4e037SJed Brown   ierr = DMSetType(*dm,DMREDUNDANT);CHKERRQ(ierr);
4728ac4e037SJed Brown   ierr = DMRedundantSetSize(*dm,rank,N);CHKERRQ(ierr);
4738ac4e037SJed Brown   ierr = DMSetUp(*dm);CHKERRQ(ierr);
4748ac4e037SJed Brown   PetscFunctionReturn(0);
4758ac4e037SJed Brown }
476