xref: /petsc/src/dm/impls/redundant/dmredundant.c (revision c688c0463f6c25494d0d243b2a43b2b3454e025c)
1b45d2f2cSJed Brown #include <petsc-private/dmimpl.h>     /*I      "petscdm.h"          I*/
28ac4e037SJed Brown #include <petscdmredundant.h>   /*I      "petscdmredundant.h" I*/
38ac4e037SJed Brown #include <petscmat.h>           /*I      "petscmat.h"         I*/
48ac4e037SJed Brown 
58ac4e037SJed Brown typedef struct  {
68ac4e037SJed Brown   PetscInt rank;                /* owner */
78ac4e037SJed Brown   PetscInt N;                   /* total number of dofs */
88ac4e037SJed Brown   PetscInt n;                   /* owned number of dofs, n=N on owner, n=0 on non-owners */
98ac4e037SJed Brown } DM_Redundant;
108ac4e037SJed Brown 
118ac4e037SJed Brown #undef __FUNCT__
12950540a4SJed Brown #define __FUNCT__ "DMCreateMatrix_Redundant"
1319fd82e9SBarry Smith static PetscErrorCode DMCreateMatrix_Redundant(DM dm,MatType mtype,Mat *J)
148ac4e037SJed Brown {
158ac4e037SJed Brown   DM_Redundant           *red = (DM_Redundant*)dm->data;
168ac4e037SJed Brown   PetscErrorCode         ierr;
178ac4e037SJed Brown   ISLocalToGlobalMapping ltog,ltogb;
18d2e2b695SJed Brown   PetscInt               i,rstart,rend,*cols;
19d2e2b695SJed Brown   PetscScalar            *vals;
208ac4e037SJed Brown 
218ac4e037SJed Brown   PetscFunctionBegin;
228ac4e037SJed Brown   ierr = MatCreate(((PetscObject)dm)->comm,J);CHKERRQ(ierr);
238ac4e037SJed Brown   ierr = MatSetSizes(*J,red->n,red->n,red->N,red->N);CHKERRQ(ierr);
248ac4e037SJed Brown   ierr = MatSetType(*J,mtype);CHKERRQ(ierr);
258ac4e037SJed Brown   ierr = MatSeqAIJSetPreallocation(*J,red->n,PETSC_NULL);CHKERRQ(ierr);
268ac4e037SJed Brown   ierr = MatSeqBAIJSetPreallocation(*J,1,red->n,PETSC_NULL);CHKERRQ(ierr);
278ac4e037SJed Brown   ierr = MatMPIAIJSetPreallocation(*J,red->n,PETSC_NULL,red->N-red->n,PETSC_NULL);CHKERRQ(ierr);
288ac4e037SJed Brown   ierr = MatMPIBAIJSetPreallocation(*J,1,red->n,PETSC_NULL,red->N-red->n,PETSC_NULL);CHKERRQ(ierr);
298ac4e037SJed Brown 
308ac4e037SJed Brown   ierr = DMGetLocalToGlobalMapping(dm,&ltog);CHKERRQ(ierr);
318ac4e037SJed Brown   ierr = DMGetLocalToGlobalMappingBlock(dm,&ltogb);CHKERRQ(ierr);
328ac4e037SJed Brown   ierr = MatSetLocalToGlobalMapping(*J,ltog,ltog);CHKERRQ(ierr);
338ac4e037SJed Brown   ierr = MatSetLocalToGlobalMappingBlock(*J,ltogb,ltogb);CHKERRQ(ierr);
34d2e2b695SJed Brown 
35d2e2b695SJed Brown   ierr = PetscMalloc2(red->N,PetscInt,&cols,red->N,PetscScalar,&vals);CHKERRQ(ierr);
36d2e2b695SJed Brown   for (i=0; i<red->N; i++) {
37d2e2b695SJed Brown     cols[i] = i;
38d2e2b695SJed Brown     vals[i] = 0.0;
39d2e2b695SJed Brown   }
40d2e2b695SJed Brown   ierr = MatGetOwnershipRange(*J,&rstart,&rend);CHKERRQ(ierr);
41d2e2b695SJed Brown   for (i=rstart; i<rend; i++) {
42d2e2b695SJed Brown     ierr = MatSetValues(*J,1,&i,red->N,cols,vals,INSERT_VALUES);CHKERRQ(ierr);
43d2e2b695SJed Brown   }
44d2e2b695SJed Brown   ierr = PetscFree2(cols,vals);CHKERRQ(ierr);
45d2e2b695SJed Brown   ierr = MatAssemblyBegin(*J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
46d2e2b695SJed Brown   ierr = MatAssemblyEnd(*J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
478ac4e037SJed Brown   PetscFunctionReturn(0);
488ac4e037SJed Brown }
498ac4e037SJed Brown 
508ac4e037SJed Brown #undef __FUNCT__
518ac4e037SJed Brown #define __FUNCT__ "DMDestroy_Redundant"
528ac4e037SJed Brown static PetscErrorCode DMDestroy_Redundant(DM dm)
538ac4e037SJed Brown {
548ac4e037SJed Brown   PetscErrorCode ierr;
558ac4e037SJed Brown 
568ac4e037SJed Brown   PetscFunctionBegin;
578ac4e037SJed Brown   ierr = PetscObjectComposeFunctionDynamic((PetscObject)dm,"DMRedundantSetSize_C","",PETSC_NULL);CHKERRQ(ierr);
588ac4e037SJed Brown   ierr = PetscObjectComposeFunctionDynamic((PetscObject)dm,"DMRedundantGetSize_C","",PETSC_NULL);CHKERRQ(ierr);
59435a35e8SMatthew G Knepley   /* This was originally freed in DMDestroy(), but that prevents reference counting of backend objects */
60435a35e8SMatthew G Knepley   ierr = PetscFree(dm->data);CHKERRQ(ierr);
618ac4e037SJed Brown   PetscFunctionReturn(0);
628ac4e037SJed Brown }
638ac4e037SJed Brown 
648ac4e037SJed Brown #undef __FUNCT__
658ac4e037SJed Brown #define __FUNCT__ "DMCreateGlobalVector_Redundant"
668ac4e037SJed Brown static PetscErrorCode DMCreateGlobalVector_Redundant(DM dm,Vec *gvec)
678ac4e037SJed Brown {
688ac4e037SJed Brown   PetscErrorCode         ierr;
698ac4e037SJed Brown   DM_Redundant           *red = (DM_Redundant*)dm->data;
708ac4e037SJed Brown   ISLocalToGlobalMapping ltog,ltogb;
718ac4e037SJed Brown 
728ac4e037SJed Brown   PetscFunctionBegin;
738ac4e037SJed Brown   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
748ac4e037SJed Brown   PetscValidPointer(gvec,2);
758ac4e037SJed Brown   *gvec = 0;
768ac4e037SJed Brown   ierr = VecCreate(((PetscObject)dm)->comm,gvec);CHKERRQ(ierr);
778ac4e037SJed Brown   ierr = VecSetSizes(*gvec,red->n,red->N);CHKERRQ(ierr);
788ac4e037SJed Brown   ierr = VecSetType(*gvec,dm->vectype);CHKERRQ(ierr);
798ac4e037SJed Brown   ierr = DMGetLocalToGlobalMapping(dm,&ltog);CHKERRQ(ierr);
808ac4e037SJed Brown   ierr = DMGetLocalToGlobalMappingBlock(dm,&ltogb);CHKERRQ(ierr);
818ac4e037SJed Brown   ierr = VecSetLocalToGlobalMapping(*gvec,ltog);CHKERRQ(ierr);
828ac4e037SJed Brown   ierr = VecSetLocalToGlobalMappingBlock(*gvec,ltog);CHKERRQ(ierr);
83*c688c046SMatthew G Knepley   ierr = VecSetDM(*gvec,dm);CHKERRQ(ierr);
848ac4e037SJed Brown   PetscFunctionReturn(0);
858ac4e037SJed Brown }
868ac4e037SJed Brown 
878ac4e037SJed Brown #undef __FUNCT__
888ac4e037SJed Brown #define __FUNCT__ "DMCreateLocalVector_Redundant"
898ac4e037SJed Brown static PetscErrorCode DMCreateLocalVector_Redundant(DM dm,Vec *lvec)
908ac4e037SJed Brown {
918ac4e037SJed Brown   PetscErrorCode ierr;
928ac4e037SJed Brown   DM_Redundant   *red = (DM_Redundant*)dm->data;
938ac4e037SJed Brown 
948ac4e037SJed Brown   PetscFunctionBegin;
958ac4e037SJed Brown   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
968ac4e037SJed Brown   PetscValidPointer(lvec,2);
978ac4e037SJed Brown   *lvec = 0;
988ac4e037SJed Brown   ierr = VecCreate(PETSC_COMM_SELF,lvec);CHKERRQ(ierr);
998ac4e037SJed Brown   ierr = VecSetSizes(*lvec,red->N,red->N);CHKERRQ(ierr);
1008ac4e037SJed Brown   ierr = VecSetType(*lvec,dm->vectype);CHKERRQ(ierr);
101*c688c046SMatthew G Knepley   ierr = VecSetDM(*lvec,dm);CHKERRQ(ierr);
1028ac4e037SJed Brown   PetscFunctionReturn(0);
1038ac4e037SJed Brown }
1048ac4e037SJed Brown 
1058ac4e037SJed Brown #undef __FUNCT__
1068ac4e037SJed Brown #define __FUNCT__ "DMLocalToGlobalBegin_Redundant"
1078ac4e037SJed Brown static PetscErrorCode DMLocalToGlobalBegin_Redundant(DM dm,Vec l,InsertMode imode,Vec g)
1088ac4e037SJed Brown {
1098ac4e037SJed Brown   PetscErrorCode    ierr;
1108ac4e037SJed Brown   DM_Redundant      *red = (DM_Redundant*)dm->data;
1118ac4e037SJed Brown   const PetscScalar *lv;
1128ac4e037SJed Brown   PetscScalar       *gv;
1138ac4e037SJed Brown   PetscMPIInt       rank;
1148ac4e037SJed Brown 
1158ac4e037SJed Brown   PetscFunctionBegin;
1168ac4e037SJed Brown   ierr = MPI_Comm_rank(((PetscObject)dm)->comm,&rank);CHKERRQ(ierr);
1178ac4e037SJed Brown   ierr = VecGetArrayRead(l,&lv);CHKERRQ(ierr);
1188ac4e037SJed Brown   ierr = VecGetArray(g,&gv);CHKERRQ(ierr);
1198ac4e037SJed Brown   switch (imode) {
1208ac4e037SJed Brown   case ADD_VALUES:
1218ac4e037SJed Brown   case MAX_VALUES:
1228ac4e037SJed Brown   {
1238ac4e037SJed Brown     void *source;
1248ac4e037SJed Brown     PetscScalar *buffer;
1258ac4e037SJed Brown     PetscInt i;
1268ac4e037SJed Brown     if (rank == red->rank) {
1278ac4e037SJed Brown #if defined (PETSC_HAVE_MPI_IN_PLACE)
1288ac4e037SJed Brown       buffer = gv;
1298ac4e037SJed Brown       source = MPI_IN_PLACE;
1308ac4e037SJed Brown #else
1318ac4e037SJed Brown       ierr = PetscMalloc(red->N*sizeof(PetscScalar),&buffer);CHKERRQ(ierr);
1328ac4e037SJed Brown       source = buffer;
1338ac4e037SJed Brown #endif
1348ac4e037SJed Brown       if (imode == ADD_VALUES) for (i=0; i<red->N; i++) buffer[i] = gv[i] + lv[i];
1358ac4e037SJed Brown #if !defined(PETSC_USE_COMPLEX)
1368ac4e037SJed Brown       if (imode == MAX_VALUES) for (i=0; i<red->N; i++) buffer[i] = PetscMax(gv[i],lv[i]);
1378ac4e037SJed Brown #endif
1388ac4e037SJed Brown     } else {
1398ac4e037SJed Brown       source = (void*)lv;
1408ac4e037SJed Brown     }
1418ac4e037SJed Brown     ierr = MPI_Reduce(source,gv,red->N,MPIU_SCALAR,(imode == ADD_VALUES)?MPI_SUM:MPI_MAX,red->rank,((PetscObject)dm)->comm);CHKERRQ(ierr);
1428ac4e037SJed Brown #if !defined(PETSC_HAVE_MPI_IN_PLACE)
14399839c1eSBarry Smith     if (rank == red->rank) {ierr = PetscFree(buffer);CHKERRQ(ierr);}
1448ac4e037SJed Brown #endif
1458ac4e037SJed Brown   } break;
1468ac4e037SJed Brown   case INSERT_VALUES:
1478ac4e037SJed Brown   ierr = PetscMemcpy(gv,lv,red->n*sizeof(PetscScalar));CHKERRQ(ierr);
1488ac4e037SJed Brown   break;
1498ac4e037SJed Brown   default: SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_SUP,"InsertMode not supported");
1508ac4e037SJed Brown   }
1518ac4e037SJed Brown   ierr = VecRestoreArrayRead(l,&lv);CHKERRQ(ierr);
1528ac4e037SJed Brown   ierr = VecRestoreArray(g,&gv);CHKERRQ(ierr);
1538ac4e037SJed Brown   PetscFunctionReturn(0);
1548ac4e037SJed Brown }
1558ac4e037SJed Brown 
1568ac4e037SJed Brown #undef __FUNCT__
1578ac4e037SJed Brown #define __FUNCT__ "DMLocalToGlobalEnd_Redundant"
1588ac4e037SJed Brown static PetscErrorCode DMLocalToGlobalEnd_Redundant(DM dm,Vec l,InsertMode imode,Vec g)
1598ac4e037SJed Brown {
1608ac4e037SJed Brown   PetscFunctionBegin;
1618ac4e037SJed Brown   PetscFunctionReturn(0);
1628ac4e037SJed Brown }
1638ac4e037SJed Brown 
1648ac4e037SJed Brown #undef __FUNCT__
1658ac4e037SJed Brown #define __FUNCT__ "DMGlobalToLocalBegin_Redundant"
1668ac4e037SJed Brown static PetscErrorCode DMGlobalToLocalBegin_Redundant(DM dm,Vec g,InsertMode imode,Vec l)
1678ac4e037SJed Brown {
1688ac4e037SJed Brown   PetscErrorCode    ierr;
1698ac4e037SJed Brown   DM_Redundant      *red = (DM_Redundant*)dm->data;
1708ac4e037SJed Brown   const PetscScalar *gv;
1718ac4e037SJed Brown   PetscScalar       *lv;
1728ac4e037SJed Brown 
1738ac4e037SJed Brown   PetscFunctionBegin;
1748ac4e037SJed Brown   ierr = VecGetArrayRead(g,&gv);CHKERRQ(ierr);
1758ac4e037SJed Brown   ierr = VecGetArray(l,&lv);CHKERRQ(ierr);
1768ac4e037SJed Brown   switch (imode) {
1778ac4e037SJed Brown   case INSERT_VALUES:
1788ac4e037SJed Brown     if (red->n) {ierr = PetscMemcpy(lv,gv,red->n*sizeof(PetscScalar));CHKERRQ(ierr);}
1798ac4e037SJed Brown     ierr = MPI_Bcast(lv,red->N,MPIU_SCALAR,red->rank,((PetscObject)dm)->comm);CHKERRQ(ierr);
1808ac4e037SJed Brown     break;
1818ac4e037SJed Brown   default: SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_SUP,"InsertMode not supported");
1828ac4e037SJed Brown   }
1838ac4e037SJed Brown   ierr = VecRestoreArrayRead(g,&gv);CHKERRQ(ierr);
1848ac4e037SJed Brown   ierr = VecRestoreArray(l,&lv);CHKERRQ(ierr);
1858ac4e037SJed Brown   PetscFunctionReturn(0);
1868ac4e037SJed Brown }
1878ac4e037SJed Brown 
1888ac4e037SJed Brown #undef __FUNCT__
1898ac4e037SJed Brown #define __FUNCT__ "DMGlobalToLocalEnd_Redundant"
1908ac4e037SJed Brown static PetscErrorCode DMGlobalToLocalEnd_Redundant(DM dm,Vec g,InsertMode imode,Vec l)
1918ac4e037SJed Brown {
1928ac4e037SJed Brown   PetscFunctionBegin;
1938ac4e037SJed Brown   PetscFunctionReturn(0);
1948ac4e037SJed Brown }
1958ac4e037SJed Brown 
1968ac4e037SJed Brown #undef __FUNCT__
1978ac4e037SJed Brown #define __FUNCT__ "DMSetUp_Redundant"
1988ac4e037SJed Brown static PetscErrorCode DMSetUp_Redundant(DM dm)
1998ac4e037SJed Brown {
2008ac4e037SJed Brown   PetscErrorCode ierr;
2018ac4e037SJed Brown   DM_Redundant   *red = (DM_Redundant*)dm->data;
2028ac4e037SJed Brown   PetscInt       i,*globals;
2038ac4e037SJed Brown 
2048ac4e037SJed Brown   PetscFunctionBegin;
2058ac4e037SJed Brown   ierr = PetscMalloc(red->N*sizeof(PetscInt),&globals);CHKERRQ(ierr);
2068ac4e037SJed Brown   for (i=0; i<red->N; i++) globals[i] = i;
2078ac4e037SJed Brown   ierr = ISLocalToGlobalMappingCreate(PETSC_COMM_SELF,red->N,globals,PETSC_OWN_POINTER,&dm->ltogmap);CHKERRQ(ierr);
2088ac4e037SJed Brown   ierr = PetscObjectReference((PetscObject)dm->ltogmap);CHKERRQ(ierr);
2098ac4e037SJed Brown   dm->ltogmapb = dm->ltogmap;
2108ac4e037SJed Brown   PetscFunctionReturn(0);
2118ac4e037SJed Brown }
2128ac4e037SJed Brown 
2138ac4e037SJed Brown #undef __FUNCT__
2148ac4e037SJed Brown #define __FUNCT__ "DMView_Redundant"
2158ac4e037SJed Brown static PetscErrorCode DMView_Redundant(DM dm,PetscViewer viewer)
2168ac4e037SJed Brown {
2178ac4e037SJed Brown   PetscErrorCode ierr;
2188ac4e037SJed Brown   DM_Redundant   *red = (DM_Redundant*)dm->data;
2198ac4e037SJed Brown   PetscBool      iascii;
2208ac4e037SJed Brown 
2218ac4e037SJed Brown   PetscFunctionBegin;
222251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
2238ac4e037SJed Brown   if (iascii) {
2248ac4e037SJed Brown     ierr = PetscViewerASCIIPrintf(viewer,"redundant: rank=%D N=%D\n",red->rank,red->N);CHKERRQ(ierr);
2258ac4e037SJed Brown   }
2268ac4e037SJed Brown   PetscFunctionReturn(0);
2278ac4e037SJed Brown }
2288ac4e037SJed Brown 
22969527181SJed Brown #undef __FUNCT__
230e727c939SJed Brown #define __FUNCT__ "DMCreateColoring_Redundant"
23119fd82e9SBarry Smith static PetscErrorCode DMCreateColoring_Redundant(DM dm,ISColoringType ctype,MatType mtype,ISColoring *coloring)
23269527181SJed Brown {
23369527181SJed Brown   DM_Redundant    *red = (DM_Redundant*)dm->data;
23469527181SJed Brown   PetscErrorCode  ierr;
23569527181SJed Brown   PetscInt        i,nloc;
23669527181SJed Brown   ISColoringValue *colors;
23769527181SJed Brown 
23869527181SJed Brown   PetscFunctionBegin;
23969527181SJed Brown   switch (ctype) {
24069527181SJed Brown   case IS_COLORING_GLOBAL:
24169527181SJed Brown     nloc = red->n;
24269527181SJed Brown     break;
24369527181SJed Brown   case IS_COLORING_GHOSTED:
24469527181SJed Brown     nloc = red->N;
24569527181SJed Brown     break;
24669527181SJed Brown   default: SETERRQ1(((PetscObject)dm)->comm,PETSC_ERR_ARG_WRONG,"Unknown ISColoringType %d",(int)ctype);
24769527181SJed Brown   }
24869527181SJed Brown   ierr = PetscMalloc(nloc*sizeof(ISColoringValue),&colors);CHKERRQ(ierr);
24969527181SJed Brown   for (i=0; i<nloc; i++) colors[i] = i;
25069527181SJed Brown   ierr = ISColoringCreate(((PetscObject)dm)->comm,red->N,nloc,colors,coloring);CHKERRQ(ierr);
25169527181SJed Brown   ierr = ISColoringSetType(*coloring,ctype);CHKERRQ(ierr);
25269527181SJed Brown   PetscFunctionReturn(0);
25369527181SJed Brown }
2548ac4e037SJed Brown 
2558ac4e037SJed Brown #undef __FUNCT__
2568ac4e037SJed Brown #define __FUNCT__ "DMRefine_Redundant"
2578ac4e037SJed Brown static PetscErrorCode DMRefine_Redundant(DM dmc,MPI_Comm comm,DM *dmf)
2588ac4e037SJed Brown {
2598ac4e037SJed Brown   PetscErrorCode ierr;
2608ac4e037SJed Brown   PetscMPIInt    flag;
2618ac4e037SJed Brown   DM_Redundant   *redc = (DM_Redundant*)dmc->data;
2628ac4e037SJed Brown 
2638ac4e037SJed Brown   PetscFunctionBegin;
2642ee06e3bSJed Brown   if (comm == MPI_COMM_NULL) comm = ((PetscObject)dmc)->comm;
2658ac4e037SJed Brown   ierr = MPI_Comm_compare(((PetscObject)dmc)->comm,comm,&flag); CHKERRQ(ierr);
2668ac4e037SJed Brown   if (flag != MPI_CONGRUENT && flag != MPI_IDENT) SETERRQ(((PetscObject)dmc)->comm,PETSC_ERR_SUP,"cannot change communicators");
2678ac4e037SJed Brown   ierr = DMRedundantCreate(comm,redc->rank,redc->N,dmf);CHKERRQ(ierr);
2688ac4e037SJed Brown   PetscFunctionReturn(0);
2698ac4e037SJed Brown }
2708ac4e037SJed Brown 
2718ac4e037SJed Brown #undef __FUNCT__
2728ac4e037SJed Brown #define __FUNCT__ "DMCoarsen_Redundant"
2738ac4e037SJed Brown static PetscErrorCode DMCoarsen_Redundant(DM dmf,MPI_Comm comm,DM *dmc)
2748ac4e037SJed Brown {
2758ac4e037SJed Brown   PetscErrorCode ierr;
2768ac4e037SJed Brown   PetscMPIInt    flag;
2778ac4e037SJed Brown   DM_Redundant   *redf = (DM_Redundant*)dmf->data;
2788ac4e037SJed Brown 
2798ac4e037SJed Brown   PetscFunctionBegin;
2802ee06e3bSJed Brown   if (comm == MPI_COMM_NULL) comm = ((PetscObject)dmf)->comm;
2818ac4e037SJed Brown   ierr = MPI_Comm_compare(((PetscObject)dmf)->comm,comm,&flag); CHKERRQ(ierr);
2828ac4e037SJed Brown   if (flag != MPI_CONGRUENT && flag != MPI_IDENT) SETERRQ(((PetscObject)dmf)->comm,PETSC_ERR_SUP,"cannot change communicators");
2838ac4e037SJed Brown   ierr = DMRedundantCreate(comm,redf->rank,redf->N,dmc);CHKERRQ(ierr);
2848ac4e037SJed Brown   PetscFunctionReturn(0);
2858ac4e037SJed Brown }
2868ac4e037SJed Brown 
2878ac4e037SJed Brown #undef __FUNCT__
288e727c939SJed Brown #define __FUNCT__ "DMCreateInterpolation_Redundant"
289e727c939SJed Brown static PetscErrorCode DMCreateInterpolation_Redundant(DM dmc,DM dmf,Mat *P,Vec *scale)
2908ac4e037SJed Brown {
2918ac4e037SJed Brown   PetscErrorCode ierr;
2928ac4e037SJed Brown   DM_Redundant   *redc = (DM_Redundant*)dmc->data;
2938ac4e037SJed Brown   DM_Redundant   *redf = (DM_Redundant*)dmf->data;
2948ac4e037SJed Brown   PetscMPIInt    flag;
2958ac4e037SJed Brown   PetscInt       i,rstart,rend;
2968ac4e037SJed Brown 
2978ac4e037SJed Brown   PetscFunctionBegin;
2988ac4e037SJed Brown   ierr = MPI_Comm_compare(((PetscObject)dmc)->comm,((PetscObject)dmf)->comm,&flag); CHKERRQ(ierr);
2998ac4e037SJed Brown   if (flag != MPI_CONGRUENT && flag != MPI_IDENT) SETERRQ(((PetscObject)dmf)->comm,PETSC_ERR_SUP,"cannot change communicators");
3008ac4e037SJed Brown   if (redc->rank != redf->rank) SETERRQ(((PetscObject)dmf)->comm,PETSC_ERR_ARG_INCOMP,"Owning rank does not match");
3018ac4e037SJed Brown   if (redc->N != redf->N) SETERRQ(((PetscObject)dmf)->comm,PETSC_ERR_ARG_INCOMP,"Global size does not match");
3028ac4e037SJed Brown   ierr = MatCreate(((PetscObject)dmc)->comm,P);CHKERRQ(ierr);
3038ac4e037SJed Brown   ierr = MatSetSizes(*P,redc->n,redc->n,redc->N,redc->N);CHKERRQ(ierr);
3048ac4e037SJed Brown   ierr = MatSetType(*P,MATAIJ);CHKERRQ(ierr);
3058ac4e037SJed Brown   ierr = MatSeqAIJSetPreallocation(*P,1,0);CHKERRQ(ierr);
3068ac4e037SJed Brown   ierr = MatMPIAIJSetPreallocation(*P,1,0,0,0);CHKERRQ(ierr);
3078ac4e037SJed Brown   ierr = MatGetOwnershipRange(*P,&rstart,&rend);CHKERRQ(ierr);
3088ac4e037SJed Brown   for (i=rstart; i<rend; i++) {ierr = MatSetValue(*P,i,i,1.0,INSERT_VALUES);CHKERRQ(ierr);}
3098ac4e037SJed Brown   ierr = MatAssemblyBegin(*P,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
3108ac4e037SJed Brown   ierr = MatAssemblyEnd(*P,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
311e727c939SJed Brown   if (scale) {ierr = DMCreateInterpolationScale(dmc,dmf,*P,scale);CHKERRQ(ierr);}
3128ac4e037SJed Brown   PetscFunctionReturn(0);
3138ac4e037SJed Brown }
3148ac4e037SJed Brown 
3158ac4e037SJed Brown #undef __FUNCT__
3168ac4e037SJed Brown #define __FUNCT__ "DMRedundantSetSize"
3178ac4e037SJed Brown /*@
3188ac4e037SJed Brown     DMRedundantSetSize - Sets the size of a densely coupled redundant object
3198ac4e037SJed Brown 
3208ac4e037SJed Brown     Collective on DM
3218ac4e037SJed Brown 
3228ac4e037SJed Brown     Input Parameter:
3238ac4e037SJed Brown +   dm - redundant DM
3248ac4e037SJed Brown .   rank - rank of process to own redundant degrees of freedom
3258ac4e037SJed Brown -   N - total number of redundant degrees of freedom
3268ac4e037SJed Brown 
3278ac4e037SJed Brown     Level: advanced
3288ac4e037SJed Brown 
3298ac4e037SJed Brown .seealso DMDestroy(), DMCreateGlobalVector(), DMRedundantCreate(), DMRedundantGetSize()
3308ac4e037SJed Brown @*/
3318ac4e037SJed Brown PetscErrorCode DMRedundantSetSize(DM dm,PetscInt rank,PetscInt N)
3328ac4e037SJed Brown {
3338ac4e037SJed Brown   PetscErrorCode ierr;
3348ac4e037SJed Brown 
3358ac4e037SJed Brown   PetscFunctionBegin;
3368ac4e037SJed Brown   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
3378ac4e037SJed Brown   PetscValidType(dm,1);
3388ac4e037SJed Brown   PetscValidLogicalCollectiveInt(dm,rank,2);
3398ac4e037SJed Brown   PetscValidLogicalCollectiveInt(dm,N,3);
3408ac4e037SJed Brown   ierr = PetscTryMethod(dm,"DMRedundantSetSize_C",(DM,PetscInt,PetscInt),(dm,rank,N));CHKERRQ(ierr);
3418ac4e037SJed Brown   PetscFunctionReturn(0);
3428ac4e037SJed Brown }
3438ac4e037SJed Brown 
3448ac4e037SJed Brown #undef __FUNCT__
3458ac4e037SJed Brown #define __FUNCT__ "DMRedundantGetSize"
3468ac4e037SJed Brown /*@
3478ac4e037SJed Brown     DMRedundantGetSize - Gets the size of a densely coupled redundant object
3488ac4e037SJed Brown 
3498ac4e037SJed Brown     Not Collective
3508ac4e037SJed Brown 
3518ac4e037SJed Brown     Input Parameter:
3528ac4e037SJed Brown +   dm - redundant DM
3538ac4e037SJed Brown 
3548ac4e037SJed Brown     Output Parameters:
3558ac4e037SJed Brown +   rank - rank of process to own redundant degrees of freedom (or PETSC_NULL)
3568ac4e037SJed Brown -   N - total number of redundant degrees of freedom (or PETSC_NULL)
3578ac4e037SJed Brown 
3588ac4e037SJed Brown     Level: advanced
3598ac4e037SJed Brown 
3608ac4e037SJed Brown .seealso DMDestroy(), DMCreateGlobalVector(), DMRedundantCreate(), DMRedundantSetSize()
3618ac4e037SJed Brown @*/
3628ac4e037SJed Brown PetscErrorCode DMRedundantGetSize(DM dm,PetscInt *rank,PetscInt *N)
3638ac4e037SJed Brown {
3648ac4e037SJed Brown   PetscErrorCode ierr;
3658ac4e037SJed Brown 
3668ac4e037SJed Brown   PetscFunctionBegin;
3678ac4e037SJed Brown   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
3688ac4e037SJed Brown   PetscValidType(dm,1);
3698ac4e037SJed Brown   ierr = PetscUseMethod(dm,"DMRedundantGetSize_C",(DM,PetscInt*,PetscInt*),(dm,rank,N));CHKERRQ(ierr);
3708ac4e037SJed Brown   PetscFunctionReturn(0);
3718ac4e037SJed Brown }
3728ac4e037SJed Brown 
3738ac4e037SJed Brown EXTERN_C_BEGIN
3748ac4e037SJed Brown #undef __FUNCT__
3758ac4e037SJed Brown #define __FUNCT__ "DMRedundantSetSize_Redundant"
3768ac4e037SJed Brown PetscErrorCode DMRedundantSetSize_Redundant(DM dm,PetscInt rank,PetscInt N)
3778ac4e037SJed Brown {
3788ac4e037SJed Brown   DM_Redundant   *red = (DM_Redundant*)dm->data;
3798ac4e037SJed Brown   PetscErrorCode ierr;
3808ac4e037SJed Brown   PetscMPIInt    myrank;
3818ac4e037SJed Brown 
3828ac4e037SJed Brown   PetscFunctionBegin;
3838ac4e037SJed Brown   ierr = MPI_Comm_rank(((PetscObject)dm)->comm,&myrank);CHKERRQ(ierr);
3848ac4e037SJed Brown   red->rank = rank;
3858ac4e037SJed Brown   red->N = N;
3868ac4e037SJed Brown   red->n = (myrank == rank) ? N : 0;
3878ac4e037SJed Brown   PetscFunctionReturn(0);
3888ac4e037SJed Brown }
3898ac4e037SJed Brown 
3908ac4e037SJed Brown #undef __FUNCT__
3918ac4e037SJed Brown #define __FUNCT__ "DMRedundantGetSize_Redundant"
3928ac4e037SJed Brown PetscErrorCode DMRedundantGetSize_Redundant(DM dm,PetscInt *rank,PetscInt *N)
3938ac4e037SJed Brown {
3948ac4e037SJed Brown   DM_Redundant   *red = (DM_Redundant*)dm->data;
3958ac4e037SJed Brown 
3968ac4e037SJed Brown   PetscFunctionBegin;
3978ac4e037SJed Brown   if (rank) *rank = red->rank;
3988ac4e037SJed Brown   if (N)    *N = red->N;
3998ac4e037SJed Brown   PetscFunctionReturn(0);
4008ac4e037SJed Brown }
4018ac4e037SJed Brown EXTERN_C_END
4028ac4e037SJed Brown 
4033efe6655SBarry Smith /*MC
4043efe6655SBarry Smith    DMREDUNDANT = "redundant" - A DM object that is used to manage data for a small set of dense globally coupled variables.
4053efe6655SBarry Smith          In the global representation of the vector the variables are all stored on a single MPI process (all the other MPI processes have
4063efe6655SBarry Smith          no variables) in the local representation all the variables are stored on ALL the MPI processes (because they are all needed for each
4073efe6655SBarry Smith          processes local computations).
4083efe6655SBarry Smith 
4093efe6655SBarry Smith          This DM is generally used inside a DMCOMPOSITE object. For example, it may be used to store continuation parameters for a bifurcation problem.
4103efe6655SBarry Smith 
4113efe6655SBarry Smith 
4123efe6655SBarry Smith   Level: intermediate
4133efe6655SBarry Smith 
4143efe6655SBarry Smith .seealso: DMType, DMCOMPOSITE, DMCreateRedundant(), DMCreate(), DMRedundantSetSize(), DMRedundantGetSize()
4153efe6655SBarry Smith M*/
4163efe6655SBarry Smith 
4178ac4e037SJed Brown EXTERN_C_BEGIN
4188ac4e037SJed Brown #undef __FUNCT__
4198ac4e037SJed Brown #define __FUNCT__ "DMCreate_Redundant"
4208ac4e037SJed Brown PetscErrorCode DMCreate_Redundant(DM dm)
4218ac4e037SJed Brown {
4228ac4e037SJed Brown   PetscErrorCode ierr;
4238ac4e037SJed Brown   DM_Redundant   *red;
4248ac4e037SJed Brown 
4258ac4e037SJed Brown   PetscFunctionBegin;
4268ac4e037SJed Brown   ierr = PetscNewLog(dm,DM_Redundant,&red);CHKERRQ(ierr);
4278ac4e037SJed Brown   dm->data = red;
4288ac4e037SJed Brown 
4298ac4e037SJed Brown   ierr = PetscObjectChangeTypeName((PetscObject)dm,DMREDUNDANT);CHKERRQ(ierr);
4308ac4e037SJed Brown   dm->ops->setup              = DMSetUp_Redundant;
4318ac4e037SJed Brown   dm->ops->view               = DMView_Redundant;
4328ac4e037SJed Brown   dm->ops->createglobalvector = DMCreateGlobalVector_Redundant;
4338ac4e037SJed Brown   dm->ops->createlocalvector  = DMCreateLocalVector_Redundant;
43425296bd5SBarry Smith   dm->ops->creatematrix       = DMCreateMatrix_Redundant;
4358ac4e037SJed Brown   dm->ops->destroy            = DMDestroy_Redundant;
4368ac4e037SJed Brown   dm->ops->globaltolocalbegin = DMGlobalToLocalBegin_Redundant;
4378ac4e037SJed Brown   dm->ops->globaltolocalend   = DMGlobalToLocalEnd_Redundant;
4388ac4e037SJed Brown   dm->ops->localtoglobalbegin = DMLocalToGlobalBegin_Redundant;
4398ac4e037SJed Brown   dm->ops->localtoglobalend   = DMLocalToGlobalEnd_Redundant;
4408ac4e037SJed Brown   dm->ops->refine             = DMRefine_Redundant;
4418ac4e037SJed Brown   dm->ops->coarsen            = DMCoarsen_Redundant;
44225296bd5SBarry Smith   dm->ops->createinterpolation= DMCreateInterpolation_Redundant;
443e727c939SJed Brown   dm->ops->getcoloring        = DMCreateColoring_Redundant;
44419fd82e9SBarry Smith   ierr = PetscStrallocpy(VECSTANDARD,(char**)&dm->vectype);CHKERRQ(ierr);
4458ac4e037SJed Brown   ierr = PetscObjectComposeFunctionDynamic((PetscObject)dm,"DMRedundantSetSize_C","DMRedundantSetSize_Redundant",DMRedundantSetSize_Redundant);CHKERRQ(ierr);
4468ac4e037SJed Brown   ierr = PetscObjectComposeFunctionDynamic((PetscObject)dm,"DMRedundantGetSize_C","DMRedundantGetSize_Redundant",DMRedundantGetSize_Redundant);CHKERRQ(ierr);
4478ac4e037SJed Brown   PetscFunctionReturn(0);
4488ac4e037SJed Brown }
4498ac4e037SJed Brown EXTERN_C_END
4508ac4e037SJed Brown 
4518ac4e037SJed Brown #undef __FUNCT__
4528ac4e037SJed Brown #define __FUNCT__ "DMRedundantCreate"
4538ac4e037SJed Brown /*@C
4548ac4e037SJed Brown     DMRedundantCreate - Creates a DM object, used to manage data for dense globally coupled variables
4558ac4e037SJed Brown 
4568ac4e037SJed Brown     Collective on MPI_Comm
4578ac4e037SJed Brown 
4588ac4e037SJed Brown     Input Parameter:
4598ac4e037SJed Brown +   comm - the processors that will share the global vector
4608ac4e037SJed Brown .   rank - rank to own the redundant values
4618ac4e037SJed Brown -   N - total number of degrees of freedom
4628ac4e037SJed Brown 
4638ac4e037SJed Brown     Output Parameters:
4648ac4e037SJed Brown .   red - the redundant DM
4658ac4e037SJed Brown 
4668ac4e037SJed Brown     Level: advanced
4678ac4e037SJed Brown 
4683efe6655SBarry Smith .seealso DMDestroy(), DMCreateGlobalVector(), DMCreateMatrix(), DMCompositeAddDM(), DMREDUNDANT, DMSetType(), DMRedundantSetSize(), DMRedundantGetSize()
4698ac4e037SJed Brown 
4708ac4e037SJed Brown @*/
4718ac4e037SJed Brown PetscErrorCode DMRedundantCreate(MPI_Comm comm,PetscInt rank,PetscInt N,DM *dm)
4728ac4e037SJed Brown {
4738ac4e037SJed Brown   PetscErrorCode ierr;
4748ac4e037SJed Brown 
4758ac4e037SJed Brown   PetscFunctionBegin;
4768ac4e037SJed Brown   PetscValidPointer(dm,2);
4778ac4e037SJed Brown   ierr = DMCreate(comm,dm);CHKERRQ(ierr);
4788ac4e037SJed Brown   ierr = DMSetType(*dm,DMREDUNDANT);CHKERRQ(ierr);
4798ac4e037SJed Brown   ierr = DMRedundantSetSize(*dm,rank,N);CHKERRQ(ierr);
4808ac4e037SJed Brown   ierr = DMSetUp(*dm);CHKERRQ(ierr);
4818ac4e037SJed Brown   PetscFunctionReturn(0);
4828ac4e037SJed Brown }
483