1*8ac4e037SJed Brown #include <private/dmimpl.h> /*I "petscdm.h" I*/ 2*8ac4e037SJed Brown #include <petscdmredundant.h> /*I "petscdmredundant.h" I*/ 3*8ac4e037SJed Brown #include <petscmat.h> /*I "petscmat.h" I*/ 4*8ac4e037SJed Brown 5*8ac4e037SJed Brown typedef struct { 6*8ac4e037SJed Brown PetscInt rank; /* owner */ 7*8ac4e037SJed Brown PetscInt N; /* total number of dofs */ 8*8ac4e037SJed Brown PetscInt n; /* owned number of dofs, n=N on owner, n=0 on non-owners */ 9*8ac4e037SJed Brown } DM_Redundant; 10*8ac4e037SJed Brown 11*8ac4e037SJed Brown #undef __FUNCT__ 12*8ac4e037SJed Brown #define __FUNCT__ "DMGetMatrix_Redundant" 13*8ac4e037SJed Brown static PetscErrorCode DMGetMatrix_Redundant(DM dm,const MatType mtype,Mat *J) 14*8ac4e037SJed Brown { 15*8ac4e037SJed Brown DM_Redundant *red = (DM_Redundant*)dm->data; 16*8ac4e037SJed Brown PetscErrorCode ierr; 17*8ac4e037SJed Brown ISLocalToGlobalMapping ltog,ltogb; 18*8ac4e037SJed Brown 19*8ac4e037SJed Brown PetscFunctionBegin; 20*8ac4e037SJed Brown ierr = MatCreate(((PetscObject)dm)->comm,J);CHKERRQ(ierr); 21*8ac4e037SJed Brown ierr = MatSetSizes(*J,red->n,red->n,red->N,red->N);CHKERRQ(ierr); 22*8ac4e037SJed Brown ierr = MatSetType(*J,mtype);CHKERRQ(ierr); 23*8ac4e037SJed Brown ierr = MatSeqAIJSetPreallocation(*J,red->n,PETSC_NULL);CHKERRQ(ierr); 24*8ac4e037SJed Brown ierr = MatSeqBAIJSetPreallocation(*J,1,red->n,PETSC_NULL);CHKERRQ(ierr); 25*8ac4e037SJed Brown ierr = MatMPIAIJSetPreallocation(*J,red->n,PETSC_NULL,red->N-red->n,PETSC_NULL);CHKERRQ(ierr); 26*8ac4e037SJed Brown ierr = MatMPIBAIJSetPreallocation(*J,1,red->n,PETSC_NULL,red->N-red->n,PETSC_NULL);CHKERRQ(ierr); 27*8ac4e037SJed Brown 28*8ac4e037SJed Brown ierr = DMGetLocalToGlobalMapping(dm,<og);CHKERRQ(ierr); 29*8ac4e037SJed Brown ierr = DMGetLocalToGlobalMappingBlock(dm,<ogb);CHKERRQ(ierr); 30*8ac4e037SJed Brown ierr = MatSetLocalToGlobalMapping(*J,ltog,ltog);CHKERRQ(ierr); 31*8ac4e037SJed Brown ierr = MatSetLocalToGlobalMappingBlock(*J,ltogb,ltogb);CHKERRQ(ierr); 32*8ac4e037SJed Brown PetscFunctionReturn(0); 33*8ac4e037SJed Brown } 34*8ac4e037SJed Brown 35*8ac4e037SJed Brown #undef __FUNCT__ 36*8ac4e037SJed Brown #define __FUNCT__ "DMDestroy_Redundant" 37*8ac4e037SJed Brown static PetscErrorCode DMDestroy_Redundant(DM dm) 38*8ac4e037SJed Brown { 39*8ac4e037SJed Brown PetscErrorCode ierr; 40*8ac4e037SJed Brown 41*8ac4e037SJed Brown PetscFunctionBegin; 42*8ac4e037SJed Brown ierr = PetscObjectComposeFunctionDynamic((PetscObject)dm,"DMRedundantSetSize_C","",PETSC_NULL);CHKERRQ(ierr); 43*8ac4e037SJed Brown ierr = PetscObjectComposeFunctionDynamic((PetscObject)dm,"DMRedundantGetSize_C","",PETSC_NULL);CHKERRQ(ierr); 44*8ac4e037SJed Brown PetscFunctionReturn(0); 45*8ac4e037SJed Brown } 46*8ac4e037SJed Brown 47*8ac4e037SJed Brown #undef __FUNCT__ 48*8ac4e037SJed Brown #define __FUNCT__ "DMCreateGlobalVector_Redundant" 49*8ac4e037SJed Brown static PetscErrorCode DMCreateGlobalVector_Redundant(DM dm,Vec *gvec) 50*8ac4e037SJed Brown { 51*8ac4e037SJed Brown PetscErrorCode ierr; 52*8ac4e037SJed Brown DM_Redundant *red = (DM_Redundant*)dm->data; 53*8ac4e037SJed Brown ISLocalToGlobalMapping ltog,ltogb; 54*8ac4e037SJed Brown 55*8ac4e037SJed Brown PetscFunctionBegin; 56*8ac4e037SJed Brown PetscValidHeaderSpecific(dm,DM_CLASSID,1); 57*8ac4e037SJed Brown PetscValidPointer(gvec,2); 58*8ac4e037SJed Brown *gvec = 0; 59*8ac4e037SJed Brown ierr = VecCreate(((PetscObject)dm)->comm,gvec);CHKERRQ(ierr); 60*8ac4e037SJed Brown ierr = VecSetSizes(*gvec,red->n,red->N);CHKERRQ(ierr); 61*8ac4e037SJed Brown ierr = VecSetType(*gvec,dm->vectype);CHKERRQ(ierr); 62*8ac4e037SJed Brown ierr = DMGetLocalToGlobalMapping(dm,<og);CHKERRQ(ierr); 63*8ac4e037SJed Brown ierr = DMGetLocalToGlobalMappingBlock(dm,<ogb);CHKERRQ(ierr); 64*8ac4e037SJed Brown ierr = VecSetLocalToGlobalMapping(*gvec,ltog);CHKERRQ(ierr); 65*8ac4e037SJed Brown ierr = VecSetLocalToGlobalMappingBlock(*gvec,ltog);CHKERRQ(ierr); 66*8ac4e037SJed Brown ierr = PetscObjectCompose((PetscObject)*gvec,"DM",(PetscObject)dm);CHKERRQ(ierr); 67*8ac4e037SJed Brown PetscFunctionReturn(0); 68*8ac4e037SJed Brown } 69*8ac4e037SJed Brown 70*8ac4e037SJed Brown #undef __FUNCT__ 71*8ac4e037SJed Brown #define __FUNCT__ "DMCreateLocalVector_Redundant" 72*8ac4e037SJed Brown static PetscErrorCode DMCreateLocalVector_Redundant(DM dm,Vec *lvec) 73*8ac4e037SJed Brown { 74*8ac4e037SJed Brown PetscErrorCode ierr; 75*8ac4e037SJed Brown DM_Redundant *red = (DM_Redundant*)dm->data; 76*8ac4e037SJed Brown 77*8ac4e037SJed Brown PetscFunctionBegin; 78*8ac4e037SJed Brown PetscValidHeaderSpecific(dm,DM_CLASSID,1); 79*8ac4e037SJed Brown PetscValidPointer(lvec,2); 80*8ac4e037SJed Brown *lvec = 0; 81*8ac4e037SJed Brown ierr = VecCreate(PETSC_COMM_SELF,lvec);CHKERRQ(ierr); 82*8ac4e037SJed Brown ierr = VecSetSizes(*lvec,red->N,red->N);CHKERRQ(ierr); 83*8ac4e037SJed Brown ierr = VecSetType(*lvec,dm->vectype);CHKERRQ(ierr); 84*8ac4e037SJed Brown ierr = PetscObjectCompose((PetscObject)*lvec,"DM",(PetscObject)dm);CHKERRQ(ierr); 85*8ac4e037SJed Brown PetscFunctionReturn(0); 86*8ac4e037SJed Brown } 87*8ac4e037SJed Brown 88*8ac4e037SJed Brown #undef __FUNCT__ 89*8ac4e037SJed Brown #define __FUNCT__ "DMLocalToGlobalBegin_Redundant" 90*8ac4e037SJed Brown static PetscErrorCode DMLocalToGlobalBegin_Redundant(DM dm,Vec l,InsertMode imode,Vec g) 91*8ac4e037SJed Brown { 92*8ac4e037SJed Brown PetscErrorCode ierr; 93*8ac4e037SJed Brown DM_Redundant *red = (DM_Redundant*)dm->data; 94*8ac4e037SJed Brown const PetscScalar *lv; 95*8ac4e037SJed Brown PetscScalar *gv; 96*8ac4e037SJed Brown PetscMPIInt rank; 97*8ac4e037SJed Brown 98*8ac4e037SJed Brown PetscFunctionBegin; 99*8ac4e037SJed Brown ierr = MPI_Comm_rank(((PetscObject)dm)->comm,&rank);CHKERRQ(ierr); 100*8ac4e037SJed Brown ierr = VecGetArrayRead(l,&lv);CHKERRQ(ierr); 101*8ac4e037SJed Brown ierr = VecGetArray(g,&gv);CHKERRQ(ierr); 102*8ac4e037SJed Brown switch (imode) { 103*8ac4e037SJed Brown case ADD_VALUES: 104*8ac4e037SJed Brown case MAX_VALUES: 105*8ac4e037SJed Brown { 106*8ac4e037SJed Brown void *source; 107*8ac4e037SJed Brown PetscScalar *buffer; 108*8ac4e037SJed Brown PetscInt i; 109*8ac4e037SJed Brown if (rank == red->rank) { 110*8ac4e037SJed Brown #if defined (PETSC_HAVE_MPI_IN_PLACE) 111*8ac4e037SJed Brown buffer = gv; 112*8ac4e037SJed Brown source = MPI_IN_PLACE; 113*8ac4e037SJed Brown #else 114*8ac4e037SJed Brown ierr = PetscMalloc(red->N*sizeof(PetscScalar),&buffer);CHKERRQ(ierr); 115*8ac4e037SJed Brown source = buffer; 116*8ac4e037SJed Brown #endif 117*8ac4e037SJed Brown if (imode == ADD_VALUES) for (i=0; i<red->N; i++) buffer[i] = gv[i] + lv[i]; 118*8ac4e037SJed Brown #if !defined(PETSC_USE_COMPLEX) 119*8ac4e037SJed Brown if (imode == MAX_VALUES) for (i=0; i<red->N; i++) buffer[i] = PetscMax(gv[i],lv[i]); 120*8ac4e037SJed Brown #endif 121*8ac4e037SJed Brown } else { 122*8ac4e037SJed Brown source = (void*)lv; 123*8ac4e037SJed Brown } 124*8ac4e037SJed Brown ierr = MPI_Reduce(source,gv,red->N,MPIU_SCALAR,(imode == ADD_VALUES)?MPI_SUM:MPI_MAX,red->rank,((PetscObject)dm)->comm);CHKERRQ(ierr); 125*8ac4e037SJed Brown #if !defined(PETSC_HAVE_MPI_IN_PLACE) 126*8ac4e037SJed Brown if (rank == mine->rank) {ierr = PetscFree(buffer);CHKERRQ(ierr);} 127*8ac4e037SJed Brown #endif 128*8ac4e037SJed Brown } break; 129*8ac4e037SJed Brown case INSERT_VALUES: 130*8ac4e037SJed Brown ierr = PetscMemcpy(gv,lv,red->n*sizeof(PetscScalar));CHKERRQ(ierr); 131*8ac4e037SJed Brown break; 132*8ac4e037SJed Brown default: SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_SUP,"InsertMode not supported"); 133*8ac4e037SJed Brown } 134*8ac4e037SJed Brown ierr = VecRestoreArrayRead(l,&lv);CHKERRQ(ierr); 135*8ac4e037SJed Brown ierr = VecRestoreArray(g,&gv);CHKERRQ(ierr); 136*8ac4e037SJed Brown PetscFunctionReturn(0); 137*8ac4e037SJed Brown } 138*8ac4e037SJed Brown 139*8ac4e037SJed Brown #undef __FUNCT__ 140*8ac4e037SJed Brown #define __FUNCT__ "DMLocalToGlobalEnd_Redundant" 141*8ac4e037SJed Brown static PetscErrorCode DMLocalToGlobalEnd_Redundant(DM dm,Vec l,InsertMode imode,Vec g) 142*8ac4e037SJed Brown { 143*8ac4e037SJed Brown PetscFunctionBegin; 144*8ac4e037SJed Brown PetscFunctionReturn(0); 145*8ac4e037SJed Brown } 146*8ac4e037SJed Brown 147*8ac4e037SJed Brown #undef __FUNCT__ 148*8ac4e037SJed Brown #define __FUNCT__ "DMGlobalToLocalBegin_Redundant" 149*8ac4e037SJed Brown static PetscErrorCode DMGlobalToLocalBegin_Redundant(DM dm,Vec g,InsertMode imode,Vec l) 150*8ac4e037SJed Brown { 151*8ac4e037SJed Brown PetscErrorCode ierr; 152*8ac4e037SJed Brown DM_Redundant *red = (DM_Redundant*)dm->data; 153*8ac4e037SJed Brown const PetscScalar *gv; 154*8ac4e037SJed Brown PetscScalar *lv; 155*8ac4e037SJed Brown 156*8ac4e037SJed Brown PetscFunctionBegin; 157*8ac4e037SJed Brown ierr = VecGetArrayRead(g,&gv);CHKERRQ(ierr); 158*8ac4e037SJed Brown ierr = VecGetArray(l,&lv);CHKERRQ(ierr); 159*8ac4e037SJed Brown switch (imode) { 160*8ac4e037SJed Brown case INSERT_VALUES: 161*8ac4e037SJed Brown if (red->n) {ierr = PetscMemcpy(lv,gv,red->n*sizeof(PetscScalar));CHKERRQ(ierr);} 162*8ac4e037SJed Brown ierr = MPI_Bcast(lv,red->N,MPIU_SCALAR,red->rank,((PetscObject)dm)->comm);CHKERRQ(ierr); 163*8ac4e037SJed Brown break; 164*8ac4e037SJed Brown default: SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_SUP,"InsertMode not supported"); 165*8ac4e037SJed Brown } 166*8ac4e037SJed Brown ierr = VecRestoreArrayRead(g,&gv);CHKERRQ(ierr); 167*8ac4e037SJed Brown ierr = VecRestoreArray(l,&lv);CHKERRQ(ierr); 168*8ac4e037SJed Brown PetscFunctionReturn(0); 169*8ac4e037SJed Brown } 170*8ac4e037SJed Brown 171*8ac4e037SJed Brown #undef __FUNCT__ 172*8ac4e037SJed Brown #define __FUNCT__ "DMGlobalToLocalEnd_Redundant" 173*8ac4e037SJed Brown static PetscErrorCode DMGlobalToLocalEnd_Redundant(DM dm,Vec g,InsertMode imode,Vec l) 174*8ac4e037SJed Brown { 175*8ac4e037SJed Brown PetscFunctionBegin; 176*8ac4e037SJed Brown PetscFunctionReturn(0); 177*8ac4e037SJed Brown } 178*8ac4e037SJed Brown 179*8ac4e037SJed Brown #undef __FUNCT__ 180*8ac4e037SJed Brown #define __FUNCT__ "DMSetUp_Redundant" 181*8ac4e037SJed Brown static PetscErrorCode DMSetUp_Redundant(DM dm) 182*8ac4e037SJed Brown { 183*8ac4e037SJed Brown PetscErrorCode ierr; 184*8ac4e037SJed Brown DM_Redundant *red = (DM_Redundant*)dm->data; 185*8ac4e037SJed Brown PetscInt i,*globals; 186*8ac4e037SJed Brown 187*8ac4e037SJed Brown PetscFunctionBegin; 188*8ac4e037SJed Brown ierr = PetscMalloc(red->N*sizeof(PetscInt),&globals);CHKERRQ(ierr); 189*8ac4e037SJed Brown for (i=0; i<red->N; i++) globals[i] = i; 190*8ac4e037SJed Brown ierr = ISLocalToGlobalMappingCreate(PETSC_COMM_SELF,red->N,globals,PETSC_OWN_POINTER,&dm->ltogmap);CHKERRQ(ierr); 191*8ac4e037SJed Brown ierr = PetscObjectReference((PetscObject)dm->ltogmap);CHKERRQ(ierr); 192*8ac4e037SJed Brown dm->ltogmapb = dm->ltogmap; 193*8ac4e037SJed Brown PetscFunctionReturn(0); 194*8ac4e037SJed Brown } 195*8ac4e037SJed Brown 196*8ac4e037SJed Brown #undef __FUNCT__ 197*8ac4e037SJed Brown #define __FUNCT__ "DMView_Redundant" 198*8ac4e037SJed Brown static PetscErrorCode DMView_Redundant(DM dm,PetscViewer viewer) 199*8ac4e037SJed Brown { 200*8ac4e037SJed Brown PetscErrorCode ierr; 201*8ac4e037SJed Brown DM_Redundant *red = (DM_Redundant*)dm->data; 202*8ac4e037SJed Brown PetscBool iascii; 203*8ac4e037SJed Brown 204*8ac4e037SJed Brown PetscFunctionBegin; 205*8ac4e037SJed Brown ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 206*8ac4e037SJed Brown if (iascii) { 207*8ac4e037SJed Brown ierr = PetscViewerASCIIPrintf(viewer,"redundant: rank=%D N=%D\n",red->rank,red->N);CHKERRQ(ierr); 208*8ac4e037SJed Brown } 209*8ac4e037SJed Brown PetscFunctionReturn(0); 210*8ac4e037SJed Brown } 211*8ac4e037SJed Brown 212*8ac4e037SJed Brown 213*8ac4e037SJed Brown #undef __FUNCT__ 214*8ac4e037SJed Brown #define __FUNCT__ "DMRefine_Redundant" 215*8ac4e037SJed Brown static PetscErrorCode DMRefine_Redundant(DM dmc,MPI_Comm comm,DM *dmf) 216*8ac4e037SJed Brown { 217*8ac4e037SJed Brown PetscErrorCode ierr; 218*8ac4e037SJed Brown PetscMPIInt flag; 219*8ac4e037SJed Brown DM_Redundant *redc = (DM_Redundant*)dmc->data; 220*8ac4e037SJed Brown 221*8ac4e037SJed Brown PetscFunctionBegin; 222*8ac4e037SJed Brown ierr = MPI_Comm_compare(((PetscObject)dmc)->comm,comm,&flag); CHKERRQ(ierr); 223*8ac4e037SJed Brown if (flag != MPI_CONGRUENT && flag != MPI_IDENT) SETERRQ(((PetscObject)dmc)->comm,PETSC_ERR_SUP,"cannot change communicators"); 224*8ac4e037SJed Brown ierr = DMRedundantCreate(comm,redc->rank,redc->N,dmf);CHKERRQ(ierr); 225*8ac4e037SJed Brown PetscFunctionReturn(0); 226*8ac4e037SJed Brown } 227*8ac4e037SJed Brown 228*8ac4e037SJed Brown #undef __FUNCT__ 229*8ac4e037SJed Brown #define __FUNCT__ "DMCoarsen_Redundant" 230*8ac4e037SJed Brown static PetscErrorCode DMCoarsen_Redundant(DM dmf,MPI_Comm comm,DM *dmc) 231*8ac4e037SJed Brown { 232*8ac4e037SJed Brown PetscErrorCode ierr; 233*8ac4e037SJed Brown PetscMPIInt flag; 234*8ac4e037SJed Brown DM_Redundant *redf = (DM_Redundant*)dmf->data; 235*8ac4e037SJed Brown 236*8ac4e037SJed Brown PetscFunctionBegin; 237*8ac4e037SJed Brown ierr = MPI_Comm_compare(((PetscObject)dmf)->comm,comm,&flag); CHKERRQ(ierr); 238*8ac4e037SJed Brown if (flag != MPI_CONGRUENT && flag != MPI_IDENT) SETERRQ(((PetscObject)dmf)->comm,PETSC_ERR_SUP,"cannot change communicators"); 239*8ac4e037SJed Brown ierr = DMRedundantCreate(comm,redf->rank,redf->N,dmc);CHKERRQ(ierr); 240*8ac4e037SJed Brown PetscFunctionReturn(0); 241*8ac4e037SJed Brown } 242*8ac4e037SJed Brown 243*8ac4e037SJed Brown #undef __FUNCT__ 244*8ac4e037SJed Brown #define __FUNCT__ "DMGetInterpolation_Redundant" 245*8ac4e037SJed Brown static PetscErrorCode DMGetInterpolation_Redundant(DM dmc,DM dmf,Mat *P,Vec *scale) 246*8ac4e037SJed Brown { 247*8ac4e037SJed Brown PetscErrorCode ierr; 248*8ac4e037SJed Brown DM_Redundant *redc = (DM_Redundant*)dmc->data; 249*8ac4e037SJed Brown DM_Redundant *redf = (DM_Redundant*)dmf->data; 250*8ac4e037SJed Brown PetscMPIInt flag; 251*8ac4e037SJed Brown PetscInt i,rstart,rend; 252*8ac4e037SJed Brown 253*8ac4e037SJed Brown PetscFunctionBegin; 254*8ac4e037SJed Brown ierr = MPI_Comm_compare(((PetscObject)dmc)->comm,((PetscObject)dmf)->comm,&flag); CHKERRQ(ierr); 255*8ac4e037SJed Brown if (flag != MPI_CONGRUENT && flag != MPI_IDENT) SETERRQ(((PetscObject)dmf)->comm,PETSC_ERR_SUP,"cannot change communicators"); 256*8ac4e037SJed Brown if (redc->rank != redf->rank) SETERRQ(((PetscObject)dmf)->comm,PETSC_ERR_ARG_INCOMP,"Owning rank does not match"); 257*8ac4e037SJed Brown if (redc->N != redf->N) SETERRQ(((PetscObject)dmf)->comm,PETSC_ERR_ARG_INCOMP,"Global size does not match"); 258*8ac4e037SJed Brown ierr = MatCreate(((PetscObject)dmc)->comm,P);CHKERRQ(ierr); 259*8ac4e037SJed Brown ierr = MatSetSizes(*P,redc->n,redc->n,redc->N,redc->N);CHKERRQ(ierr); 260*8ac4e037SJed Brown ierr = MatSetType(*P,MATAIJ);CHKERRQ(ierr); 261*8ac4e037SJed Brown ierr = MatSeqAIJSetPreallocation(*P,1,0);CHKERRQ(ierr); 262*8ac4e037SJed Brown ierr = MatMPIAIJSetPreallocation(*P,1,0,0,0);CHKERRQ(ierr); 263*8ac4e037SJed Brown ierr = MatGetOwnershipRange(*P,&rstart,&rend);CHKERRQ(ierr); 264*8ac4e037SJed Brown for (i=rstart; i<rend; i++) {ierr = MatSetValue(*P,i,i,1.0,INSERT_VALUES);CHKERRQ(ierr);} 265*8ac4e037SJed Brown ierr = MatAssemblyBegin(*P,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 266*8ac4e037SJed Brown ierr = MatAssemblyEnd(*P,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 267*8ac4e037SJed Brown if (scale) {ierr = DMGetInterpolationScale(dmc,dmf,*P,scale);CHKERRQ(ierr);} 268*8ac4e037SJed Brown PetscFunctionReturn(0); 269*8ac4e037SJed Brown } 270*8ac4e037SJed Brown 271*8ac4e037SJed Brown #undef __FUNCT__ 272*8ac4e037SJed Brown #define __FUNCT__ "DMRedundantSetSize" 273*8ac4e037SJed Brown /*@ 274*8ac4e037SJed Brown DMRedundantSetSize - Sets the size of a densely coupled redundant object 275*8ac4e037SJed Brown 276*8ac4e037SJed Brown Collective on DM 277*8ac4e037SJed Brown 278*8ac4e037SJed Brown Input Parameter: 279*8ac4e037SJed Brown + dm - redundant DM 280*8ac4e037SJed Brown . rank - rank of process to own redundant degrees of freedom 281*8ac4e037SJed Brown - N - total number of redundant degrees of freedom 282*8ac4e037SJed Brown 283*8ac4e037SJed Brown Level: advanced 284*8ac4e037SJed Brown 285*8ac4e037SJed Brown .seealso DMDestroy(), DMCreateGlobalVector(), DMRedundantCreate(), DMRedundantGetSize() 286*8ac4e037SJed Brown @*/ 287*8ac4e037SJed Brown PetscErrorCode DMRedundantSetSize(DM dm,PetscInt rank,PetscInt N) 288*8ac4e037SJed Brown { 289*8ac4e037SJed Brown PetscErrorCode ierr; 290*8ac4e037SJed Brown 291*8ac4e037SJed Brown PetscFunctionBegin; 292*8ac4e037SJed Brown PetscValidHeaderSpecific(dm,DM_CLASSID,1); 293*8ac4e037SJed Brown PetscValidType(dm,1); 294*8ac4e037SJed Brown PetscValidLogicalCollectiveInt(dm,rank,2); 295*8ac4e037SJed Brown PetscValidLogicalCollectiveInt(dm,N,3); 296*8ac4e037SJed Brown ierr = PetscTryMethod(dm,"DMRedundantSetSize_C",(DM,PetscInt,PetscInt),(dm,rank,N));CHKERRQ(ierr); 297*8ac4e037SJed Brown PetscFunctionReturn(0); 298*8ac4e037SJed Brown } 299*8ac4e037SJed Brown 300*8ac4e037SJed Brown #undef __FUNCT__ 301*8ac4e037SJed Brown #define __FUNCT__ "DMRedundantGetSize" 302*8ac4e037SJed Brown /*@ 303*8ac4e037SJed Brown DMRedundantGetSize - Gets the size of a densely coupled redundant object 304*8ac4e037SJed Brown 305*8ac4e037SJed Brown Not Collective 306*8ac4e037SJed Brown 307*8ac4e037SJed Brown Input Parameter: 308*8ac4e037SJed Brown + dm - redundant DM 309*8ac4e037SJed Brown 310*8ac4e037SJed Brown Output Parameters: 311*8ac4e037SJed Brown + rank - rank of process to own redundant degrees of freedom (or PETSC_NULL) 312*8ac4e037SJed Brown - N - total number of redundant degrees of freedom (or PETSC_NULL) 313*8ac4e037SJed Brown 314*8ac4e037SJed Brown Level: advanced 315*8ac4e037SJed Brown 316*8ac4e037SJed Brown .seealso DMDestroy(), DMCreateGlobalVector(), DMRedundantCreate(), DMRedundantSetSize() 317*8ac4e037SJed Brown @*/ 318*8ac4e037SJed Brown PetscErrorCode DMRedundantGetSize(DM dm,PetscInt *rank,PetscInt *N) 319*8ac4e037SJed Brown { 320*8ac4e037SJed Brown PetscErrorCode ierr; 321*8ac4e037SJed Brown 322*8ac4e037SJed Brown PetscFunctionBegin; 323*8ac4e037SJed Brown PetscValidHeaderSpecific(dm,DM_CLASSID,1); 324*8ac4e037SJed Brown PetscValidType(dm,1); 325*8ac4e037SJed Brown ierr = PetscUseMethod(dm,"DMRedundantGetSize_C",(DM,PetscInt*,PetscInt*),(dm,rank,N));CHKERRQ(ierr); 326*8ac4e037SJed Brown PetscFunctionReturn(0); 327*8ac4e037SJed Brown } 328*8ac4e037SJed Brown 329*8ac4e037SJed Brown EXTERN_C_BEGIN 330*8ac4e037SJed Brown #undef __FUNCT__ 331*8ac4e037SJed Brown #define __FUNCT__ "DMRedundantSetSize_Redundant" 332*8ac4e037SJed Brown PetscErrorCode DMRedundantSetSize_Redundant(DM dm,PetscInt rank,PetscInt N) 333*8ac4e037SJed Brown { 334*8ac4e037SJed Brown DM_Redundant *red = (DM_Redundant*)dm->data; 335*8ac4e037SJed Brown PetscErrorCode ierr; 336*8ac4e037SJed Brown PetscMPIInt myrank; 337*8ac4e037SJed Brown 338*8ac4e037SJed Brown PetscFunctionBegin; 339*8ac4e037SJed Brown ierr = MPI_Comm_rank(((PetscObject)dm)->comm,&myrank);CHKERRQ(ierr); 340*8ac4e037SJed Brown red->rank = rank; 341*8ac4e037SJed Brown red->N = N; 342*8ac4e037SJed Brown red->n = (myrank == rank) ? N : 0; 343*8ac4e037SJed Brown PetscFunctionReturn(0); 344*8ac4e037SJed Brown } 345*8ac4e037SJed Brown 346*8ac4e037SJed Brown #undef __FUNCT__ 347*8ac4e037SJed Brown #define __FUNCT__ "DMRedundantGetSize_Redundant" 348*8ac4e037SJed Brown PetscErrorCode DMRedundantGetSize_Redundant(DM dm,PetscInt *rank,PetscInt *N) 349*8ac4e037SJed Brown { 350*8ac4e037SJed Brown DM_Redundant *red = (DM_Redundant*)dm->data; 351*8ac4e037SJed Brown 352*8ac4e037SJed Brown PetscFunctionBegin; 353*8ac4e037SJed Brown if (rank) *rank = red->rank; 354*8ac4e037SJed Brown if (N) *N = red->N; 355*8ac4e037SJed Brown PetscFunctionReturn(0); 356*8ac4e037SJed Brown } 357*8ac4e037SJed Brown EXTERN_C_END 358*8ac4e037SJed Brown 359*8ac4e037SJed Brown EXTERN_C_BEGIN 360*8ac4e037SJed Brown #undef __FUNCT__ 361*8ac4e037SJed Brown #define __FUNCT__ "DMCreate_Redundant" 362*8ac4e037SJed Brown PetscErrorCode DMCreate_Redundant(DM dm) 363*8ac4e037SJed Brown { 364*8ac4e037SJed Brown PetscErrorCode ierr; 365*8ac4e037SJed Brown DM_Redundant *red; 366*8ac4e037SJed Brown 367*8ac4e037SJed Brown PetscFunctionBegin; 368*8ac4e037SJed Brown ierr = PetscNewLog(dm,DM_Redundant,&red);CHKERRQ(ierr); 369*8ac4e037SJed Brown dm->data = red; 370*8ac4e037SJed Brown 371*8ac4e037SJed Brown ierr = PetscObjectChangeTypeName((PetscObject)dm,DMREDUNDANT);CHKERRQ(ierr); 372*8ac4e037SJed Brown dm->ops->setup = DMSetUp_Redundant; 373*8ac4e037SJed Brown dm->ops->view = DMView_Redundant; 374*8ac4e037SJed Brown dm->ops->createglobalvector = DMCreateGlobalVector_Redundant; 375*8ac4e037SJed Brown dm->ops->createlocalvector = DMCreateLocalVector_Redundant; 376*8ac4e037SJed Brown dm->ops->getmatrix = DMGetMatrix_Redundant; 377*8ac4e037SJed Brown dm->ops->destroy = DMDestroy_Redundant; 378*8ac4e037SJed Brown dm->ops->globaltolocalbegin = DMGlobalToLocalBegin_Redundant; 379*8ac4e037SJed Brown dm->ops->globaltolocalend = DMGlobalToLocalEnd_Redundant; 380*8ac4e037SJed Brown dm->ops->localtoglobalbegin = DMLocalToGlobalBegin_Redundant; 381*8ac4e037SJed Brown dm->ops->localtoglobalend = DMLocalToGlobalEnd_Redundant; 382*8ac4e037SJed Brown dm->ops->refine = DMRefine_Redundant; 383*8ac4e037SJed Brown dm->ops->coarsen = DMCoarsen_Redundant; 384*8ac4e037SJed Brown dm->ops->getinterpolation = DMGetInterpolation_Redundant; 385*8ac4e037SJed Brown ierr = PetscStrallocpy(VECSTANDARD,&dm->vectype);CHKERRQ(ierr); 386*8ac4e037SJed Brown ierr = PetscObjectComposeFunctionDynamic((PetscObject)dm,"DMRedundantSetSize_C","DMRedundantSetSize_Redundant",DMRedundantSetSize_Redundant);CHKERRQ(ierr); 387*8ac4e037SJed Brown ierr = PetscObjectComposeFunctionDynamic((PetscObject)dm,"DMRedundantGetSize_C","DMRedundantGetSize_Redundant",DMRedundantGetSize_Redundant);CHKERRQ(ierr); 388*8ac4e037SJed Brown PetscFunctionReturn(0); 389*8ac4e037SJed Brown } 390*8ac4e037SJed Brown EXTERN_C_END 391*8ac4e037SJed Brown 392*8ac4e037SJed Brown #undef __FUNCT__ 393*8ac4e037SJed Brown #define __FUNCT__ "DMRedundantCreate" 394*8ac4e037SJed Brown /*@C 395*8ac4e037SJed Brown DMRedundantCreate - Creates a DM object, used to manage data for dense globally coupled variables 396*8ac4e037SJed Brown 397*8ac4e037SJed Brown Collective on MPI_Comm 398*8ac4e037SJed Brown 399*8ac4e037SJed Brown Input Parameter: 400*8ac4e037SJed Brown + comm - the processors that will share the global vector 401*8ac4e037SJed Brown . rank - rank to own the redundant values 402*8ac4e037SJed Brown - N - total number of degrees of freedom 403*8ac4e037SJed Brown 404*8ac4e037SJed Brown Output Parameters: 405*8ac4e037SJed Brown . red - the redundant DM 406*8ac4e037SJed Brown 407*8ac4e037SJed Brown Level: advanced 408*8ac4e037SJed Brown 409*8ac4e037SJed Brown .seealso DMDestroy(), DMCreateGlobalVector(), DMGetMatrix(), DMCompositeAddDM() 410*8ac4e037SJed Brown 411*8ac4e037SJed Brown @*/ 412*8ac4e037SJed Brown PetscErrorCode DMRedundantCreate(MPI_Comm comm,PetscInt rank,PetscInt N,DM *dm) 413*8ac4e037SJed Brown { 414*8ac4e037SJed Brown PetscErrorCode ierr; 415*8ac4e037SJed Brown 416*8ac4e037SJed Brown PetscFunctionBegin; 417*8ac4e037SJed Brown PetscValidPointer(dm,2); 418*8ac4e037SJed Brown ierr = DMCreate(comm,dm);CHKERRQ(ierr); 419*8ac4e037SJed Brown ierr = DMSetType(*dm,DMREDUNDANT);CHKERRQ(ierr); 420*8ac4e037SJed Brown ierr = DMRedundantSetSize(*dm,rank,N);CHKERRQ(ierr); 421*8ac4e037SJed Brown ierr = DMSetUp(*dm);CHKERRQ(ierr); 422*8ac4e037SJed Brown PetscFunctionReturn(0); 423*8ac4e037SJed Brown } 424