1 #include <petsc-private/dmimpl.h> 2 #include <petscdmredundant.h> /*I "petscdmredundant.h" I*/ 3 4 typedef struct { 5 PetscInt rank; /* owner */ 6 PetscInt N; /* total number of dofs */ 7 PetscInt n; /* owned number of dofs, n=N on owner, n=0 on non-owners */ 8 } DM_Redundant; 9 10 #undef __FUNCT__ 11 #define __FUNCT__ "DMCreateMatrix_Redundant" 12 static PetscErrorCode DMCreateMatrix_Redundant(DM dm,Mat *J) 13 { 14 DM_Redundant *red = (DM_Redundant*)dm->data; 15 PetscErrorCode ierr; 16 ISLocalToGlobalMapping ltog; 17 PetscInt i,rstart,rend,*cols; 18 PetscScalar *vals; 19 20 PetscFunctionBegin; 21 ierr = MatCreate(PetscObjectComm((PetscObject)dm),J);CHKERRQ(ierr); 22 ierr = MatSetSizes(*J,red->n,red->n,red->N,red->N);CHKERRQ(ierr); 23 ierr = MatSetType(*J,dm->mattype);CHKERRQ(ierr); 24 ierr = MatSeqAIJSetPreallocation(*J,red->n,NULL);CHKERRQ(ierr); 25 ierr = MatSeqBAIJSetPreallocation(*J,1,red->n,NULL);CHKERRQ(ierr); 26 ierr = MatMPIAIJSetPreallocation(*J,red->n,NULL,red->N-red->n,NULL);CHKERRQ(ierr); 27 ierr = MatMPIBAIJSetPreallocation(*J,1,red->n,NULL,red->N-red->n,NULL);CHKERRQ(ierr); 28 29 ierr = DMGetLocalToGlobalMapping(dm,<og);CHKERRQ(ierr); 30 ierr = MatSetLocalToGlobalMapping(*J,ltog,ltog);CHKERRQ(ierr); 31 32 ierr = PetscMalloc2(red->N,&cols,red->N,&vals);CHKERRQ(ierr); 33 for (i=0; i<red->N; i++) { 34 cols[i] = i; 35 vals[i] = 0.0; 36 } 37 ierr = MatGetOwnershipRange(*J,&rstart,&rend);CHKERRQ(ierr); 38 for (i=rstart; i<rend; i++) { 39 ierr = MatSetValues(*J,1,&i,red->N,cols,vals,INSERT_VALUES);CHKERRQ(ierr); 40 } 41 ierr = PetscFree2(cols,vals);CHKERRQ(ierr); 42 ierr = MatAssemblyBegin(*J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 43 ierr = MatAssemblyEnd(*J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 44 PetscFunctionReturn(0); 45 } 46 47 #undef __FUNCT__ 48 #define __FUNCT__ "DMDestroy_Redundant" 49 static PetscErrorCode DMDestroy_Redundant(DM dm) 50 { 51 PetscErrorCode ierr; 52 53 PetscFunctionBegin; 54 ierr = PetscObjectComposeFunction((PetscObject)dm,"DMRedundantSetSize_C",NULL);CHKERRQ(ierr); 55 ierr = PetscObjectComposeFunction((PetscObject)dm,"DMRedundantGetSize_C",NULL);CHKERRQ(ierr); 56 /* This was originally freed in DMDestroy(), but that prevents reference counting of backend objects */ 57 ierr = PetscFree(dm->data);CHKERRQ(ierr); 58 PetscFunctionReturn(0); 59 } 60 61 #undef __FUNCT__ 62 #define __FUNCT__ "DMCreateGlobalVector_Redundant" 63 static PetscErrorCode DMCreateGlobalVector_Redundant(DM dm,Vec *gvec) 64 { 65 PetscErrorCode ierr; 66 DM_Redundant *red = (DM_Redundant*)dm->data; 67 ISLocalToGlobalMapping ltog; 68 69 PetscFunctionBegin; 70 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 71 PetscValidPointer(gvec,2); 72 *gvec = 0; 73 ierr = VecCreate(PetscObjectComm((PetscObject)dm),gvec);CHKERRQ(ierr); 74 ierr = VecSetSizes(*gvec,red->n,red->N);CHKERRQ(ierr); 75 ierr = VecSetType(*gvec,dm->vectype);CHKERRQ(ierr); 76 ierr = DMGetLocalToGlobalMapping(dm,<og);CHKERRQ(ierr); 77 ierr = VecSetLocalToGlobalMapping(*gvec,ltog);CHKERRQ(ierr); 78 ierr = VecSetDM(*gvec,dm);CHKERRQ(ierr); 79 PetscFunctionReturn(0); 80 } 81 82 #undef __FUNCT__ 83 #define __FUNCT__ "DMCreateLocalVector_Redundant" 84 static PetscErrorCode DMCreateLocalVector_Redundant(DM dm,Vec *lvec) 85 { 86 PetscErrorCode ierr; 87 DM_Redundant *red = (DM_Redundant*)dm->data; 88 89 PetscFunctionBegin; 90 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 91 PetscValidPointer(lvec,2); 92 *lvec = 0; 93 ierr = VecCreate(PETSC_COMM_SELF,lvec);CHKERRQ(ierr); 94 ierr = VecSetSizes(*lvec,red->N,red->N);CHKERRQ(ierr); 95 ierr = VecSetType(*lvec,dm->vectype);CHKERRQ(ierr); 96 ierr = VecSetDM(*lvec,dm);CHKERRQ(ierr); 97 PetscFunctionReturn(0); 98 } 99 100 #undef __FUNCT__ 101 #define __FUNCT__ "DMLocalToGlobalBegin_Redundant" 102 static PetscErrorCode DMLocalToGlobalBegin_Redundant(DM dm,Vec l,InsertMode imode,Vec g) 103 { 104 PetscErrorCode ierr; 105 DM_Redundant *red = (DM_Redundant*)dm->data; 106 const PetscScalar *lv; 107 PetscScalar *gv; 108 PetscMPIInt rank; 109 110 PetscFunctionBegin; 111 ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)dm),&rank);CHKERRQ(ierr); 112 ierr = VecGetArrayRead(l,&lv);CHKERRQ(ierr); 113 ierr = VecGetArray(g,&gv);CHKERRQ(ierr); 114 switch (imode) { 115 case ADD_VALUES: 116 case MAX_VALUES: 117 { 118 void *source; 119 PetscScalar *buffer; 120 PetscInt i; 121 if (rank == red->rank) { 122 #if defined(PETSC_HAVE_MPI_IN_PLACE) 123 buffer = gv; 124 source = MPI_IN_PLACE; 125 #else 126 ierr = PetscMalloc1(red->N,&buffer);CHKERRQ(ierr); 127 source = buffer; 128 #endif 129 if (imode == ADD_VALUES) for (i=0; i<red->N; i++) buffer[i] = gv[i] + lv[i]; 130 #if !defined(PETSC_USE_COMPLEX) 131 if (imode == MAX_VALUES) for (i=0; i<red->N; i++) buffer[i] = PetscMax(gv[i],lv[i]); 132 #endif 133 } else source = (void*)lv; 134 ierr = MPI_Reduce(source,gv,red->N,MPIU_SCALAR,(imode == ADD_VALUES) ? MPIU_SUM : MPIU_MAX,red->rank,PetscObjectComm((PetscObject)dm));CHKERRQ(ierr); 135 #if !defined(PETSC_HAVE_MPI_IN_PLACE) 136 if (rank == red->rank) {ierr = PetscFree(buffer);CHKERRQ(ierr);} 137 #endif 138 } break; 139 case INSERT_VALUES: 140 ierr = PetscMemcpy(gv,lv,red->n*sizeof(PetscScalar));CHKERRQ(ierr); 141 break; 142 default: SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"InsertMode not supported"); 143 } 144 ierr = VecRestoreArrayRead(l,&lv);CHKERRQ(ierr); 145 ierr = VecRestoreArray(g,&gv);CHKERRQ(ierr); 146 PetscFunctionReturn(0); 147 } 148 149 #undef __FUNCT__ 150 #define __FUNCT__ "DMLocalToGlobalEnd_Redundant" 151 static PetscErrorCode DMLocalToGlobalEnd_Redundant(DM dm,Vec l,InsertMode imode,Vec g) 152 { 153 PetscFunctionBegin; 154 PetscFunctionReturn(0); 155 } 156 157 #undef __FUNCT__ 158 #define __FUNCT__ "DMGlobalToLocalBegin_Redundant" 159 static PetscErrorCode DMGlobalToLocalBegin_Redundant(DM dm,Vec g,InsertMode imode,Vec l) 160 { 161 PetscErrorCode ierr; 162 DM_Redundant *red = (DM_Redundant*)dm->data; 163 const PetscScalar *gv; 164 PetscScalar *lv; 165 166 PetscFunctionBegin; 167 ierr = VecGetArrayRead(g,&gv);CHKERRQ(ierr); 168 ierr = VecGetArray(l,&lv);CHKERRQ(ierr); 169 switch (imode) { 170 case INSERT_VALUES: 171 if (red->n) {ierr = PetscMemcpy(lv,gv,red->n*sizeof(PetscScalar));CHKERRQ(ierr);} 172 ierr = MPI_Bcast(lv,red->N,MPIU_SCALAR,red->rank,PetscObjectComm((PetscObject)dm));CHKERRQ(ierr); 173 break; 174 default: SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"InsertMode not supported"); 175 } 176 ierr = VecRestoreArrayRead(g,&gv);CHKERRQ(ierr); 177 ierr = VecRestoreArray(l,&lv);CHKERRQ(ierr); 178 PetscFunctionReturn(0); 179 } 180 181 #undef __FUNCT__ 182 #define __FUNCT__ "DMGlobalToLocalEnd_Redundant" 183 static PetscErrorCode DMGlobalToLocalEnd_Redundant(DM dm,Vec g,InsertMode imode,Vec l) 184 { 185 PetscFunctionBegin; 186 PetscFunctionReturn(0); 187 } 188 189 #undef __FUNCT__ 190 #define __FUNCT__ "DMSetUp_Redundant" 191 static PetscErrorCode DMSetUp_Redundant(DM dm) 192 { 193 PetscErrorCode ierr; 194 DM_Redundant *red = (DM_Redundant*)dm->data; 195 PetscInt i,*globals; 196 197 PetscFunctionBegin; 198 ierr = PetscMalloc1(red->N,&globals);CHKERRQ(ierr); 199 for (i=0; i<red->N; i++) globals[i] = i; 200 ierr = ISLocalToGlobalMappingCreate(PETSC_COMM_SELF,1,red->N,globals,PETSC_OWN_POINTER,&dm->ltogmap);CHKERRQ(ierr); 201 ierr = PetscObjectReference((PetscObject)dm->ltogmap);CHKERRQ(ierr); 202 PetscFunctionReturn(0); 203 } 204 205 #undef __FUNCT__ 206 #define __FUNCT__ "DMView_Redundant" 207 static PetscErrorCode DMView_Redundant(DM dm,PetscViewer viewer) 208 { 209 PetscErrorCode ierr; 210 DM_Redundant *red = (DM_Redundant*)dm->data; 211 PetscBool iascii; 212 213 PetscFunctionBegin; 214 ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 215 if (iascii) { 216 ierr = PetscViewerASCIIPrintf(viewer,"redundant: rank=%D N=%D\n",red->rank,red->N);CHKERRQ(ierr); 217 } 218 PetscFunctionReturn(0); 219 } 220 221 #undef __FUNCT__ 222 #define __FUNCT__ "DMCreateColoring_Redundant" 223 static PetscErrorCode DMCreateColoring_Redundant(DM dm,ISColoringType ctype,ISColoring *coloring) 224 { 225 DM_Redundant *red = (DM_Redundant*)dm->data; 226 PetscErrorCode ierr; 227 PetscInt i,nloc; 228 ISColoringValue *colors; 229 230 PetscFunctionBegin; 231 switch (ctype) { 232 case IS_COLORING_GLOBAL: 233 nloc = red->n; 234 break; 235 case IS_COLORING_GHOSTED: 236 nloc = red->N; 237 break; 238 default: SETERRQ1(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONG,"Unknown ISColoringType %d",(int)ctype); 239 } 240 ierr = PetscMalloc1(nloc,&colors);CHKERRQ(ierr); 241 for (i=0; i<nloc; i++) colors[i] = i; 242 ierr = ISColoringCreate(PetscObjectComm((PetscObject)dm),red->N,nloc,colors,coloring);CHKERRQ(ierr); 243 ierr = ISColoringSetType(*coloring,ctype);CHKERRQ(ierr); 244 PetscFunctionReturn(0); 245 } 246 247 #undef __FUNCT__ 248 #define __FUNCT__ "DMRefine_Redundant" 249 static PetscErrorCode DMRefine_Redundant(DM dmc,MPI_Comm comm,DM *dmf) 250 { 251 PetscErrorCode ierr; 252 PetscMPIInt flag; 253 DM_Redundant *redc = (DM_Redundant*)dmc->data; 254 255 PetscFunctionBegin; 256 if (comm == MPI_COMM_NULL) { 257 ierr = PetscObjectGetComm((PetscObject)dmc,&comm);CHKERRQ(ierr); 258 } 259 ierr = MPI_Comm_compare(PetscObjectComm((PetscObject)dmc),comm,&flag);CHKERRQ(ierr); 260 if (flag != MPI_CONGRUENT && flag != MPI_IDENT) SETERRQ(PetscObjectComm((PetscObject)dmc),PETSC_ERR_SUP,"cannot change communicators"); 261 ierr = DMRedundantCreate(comm,redc->rank,redc->N,dmf);CHKERRQ(ierr); 262 PetscFunctionReturn(0); 263 } 264 265 #undef __FUNCT__ 266 #define __FUNCT__ "DMCoarsen_Redundant" 267 static PetscErrorCode DMCoarsen_Redundant(DM dmf,MPI_Comm comm,DM *dmc) 268 { 269 PetscErrorCode ierr; 270 PetscMPIInt flag; 271 DM_Redundant *redf = (DM_Redundant*)dmf->data; 272 273 PetscFunctionBegin; 274 if (comm == MPI_COMM_NULL) { 275 ierr = PetscObjectGetComm((PetscObject)dmf,&comm);CHKERRQ(ierr); 276 } 277 ierr = MPI_Comm_compare(PetscObjectComm((PetscObject)dmf),comm,&flag);CHKERRQ(ierr); 278 if (flag != MPI_CONGRUENT && flag != MPI_IDENT) SETERRQ(PetscObjectComm((PetscObject)dmf),PETSC_ERR_SUP,"cannot change communicators"); 279 ierr = DMRedundantCreate(comm,redf->rank,redf->N,dmc);CHKERRQ(ierr); 280 PetscFunctionReturn(0); 281 } 282 283 #undef __FUNCT__ 284 #define __FUNCT__ "DMCreateInterpolation_Redundant" 285 static PetscErrorCode DMCreateInterpolation_Redundant(DM dmc,DM dmf,Mat *P,Vec *scale) 286 { 287 PetscErrorCode ierr; 288 DM_Redundant *redc = (DM_Redundant*)dmc->data; 289 DM_Redundant *redf = (DM_Redundant*)dmf->data; 290 PetscMPIInt flag; 291 PetscInt i,rstart,rend; 292 293 PetscFunctionBegin; 294 ierr = MPI_Comm_compare(PetscObjectComm((PetscObject)dmc),PetscObjectComm((PetscObject)dmf),&flag);CHKERRQ(ierr); 295 if (flag != MPI_CONGRUENT && flag != MPI_IDENT) SETERRQ(PetscObjectComm((PetscObject)dmf),PETSC_ERR_SUP,"cannot change communicators"); 296 if (redc->rank != redf->rank) SETERRQ(PetscObjectComm((PetscObject)dmf),PETSC_ERR_ARG_INCOMP,"Owning rank does not match"); 297 if (redc->N != redf->N) SETERRQ(PetscObjectComm((PetscObject)dmf),PETSC_ERR_ARG_INCOMP,"Global size does not match"); 298 ierr = MatCreate(PetscObjectComm((PetscObject)dmc),P);CHKERRQ(ierr); 299 ierr = MatSetSizes(*P,redc->n,redc->n,redc->N,redc->N);CHKERRQ(ierr); 300 ierr = MatSetType(*P,MATAIJ);CHKERRQ(ierr); 301 ierr = MatSeqAIJSetPreallocation(*P,1,0);CHKERRQ(ierr); 302 ierr = MatMPIAIJSetPreallocation(*P,1,0,0,0);CHKERRQ(ierr); 303 ierr = MatGetOwnershipRange(*P,&rstart,&rend);CHKERRQ(ierr); 304 for (i=rstart; i<rend; i++) {ierr = MatSetValue(*P,i,i,1.0,INSERT_VALUES);CHKERRQ(ierr);} 305 ierr = MatAssemblyBegin(*P,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 306 ierr = MatAssemblyEnd(*P,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 307 if (scale) {ierr = DMCreateInterpolationScale(dmc,dmf,*P,scale);CHKERRQ(ierr);} 308 PetscFunctionReturn(0); 309 } 310 311 #undef __FUNCT__ 312 #define __FUNCT__ "DMRedundantSetSize" 313 /*@ 314 DMRedundantSetSize - Sets the size of a densely coupled redundant object 315 316 Collective on DM 317 318 Input Parameter: 319 + dm - redundant DM 320 . rank - rank of process to own redundant degrees of freedom 321 - N - total number of redundant degrees of freedom 322 323 Level: advanced 324 325 .seealso DMDestroy(), DMCreateGlobalVector(), DMRedundantCreate(), DMRedundantGetSize() 326 @*/ 327 PetscErrorCode DMRedundantSetSize(DM dm,PetscInt rank,PetscInt N) 328 { 329 PetscErrorCode ierr; 330 331 PetscFunctionBegin; 332 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 333 PetscValidType(dm,1); 334 PetscValidLogicalCollectiveInt(dm,rank,2); 335 PetscValidLogicalCollectiveInt(dm,N,3); 336 ierr = PetscTryMethod(dm,"DMRedundantSetSize_C",(DM,PetscInt,PetscInt),(dm,rank,N));CHKERRQ(ierr); 337 PetscFunctionReturn(0); 338 } 339 340 #undef __FUNCT__ 341 #define __FUNCT__ "DMRedundantGetSize" 342 /*@ 343 DMRedundantGetSize - Gets the size of a densely coupled redundant object 344 345 Not Collective 346 347 Input Parameter: 348 + dm - redundant DM 349 350 Output Parameters: 351 + rank - rank of process to own redundant degrees of freedom (or NULL) 352 - N - total number of redundant degrees of freedom (or NULL) 353 354 Level: advanced 355 356 .seealso DMDestroy(), DMCreateGlobalVector(), DMRedundantCreate(), DMRedundantSetSize() 357 @*/ 358 PetscErrorCode DMRedundantGetSize(DM dm,PetscInt *rank,PetscInt *N) 359 { 360 PetscErrorCode ierr; 361 362 PetscFunctionBegin; 363 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 364 PetscValidType(dm,1); 365 ierr = PetscUseMethod(dm,"DMRedundantGetSize_C",(DM,PetscInt*,PetscInt*),(dm,rank,N));CHKERRQ(ierr); 366 PetscFunctionReturn(0); 367 } 368 369 #undef __FUNCT__ 370 #define __FUNCT__ "DMRedundantSetSize_Redundant" 371 static PetscErrorCode DMRedundantSetSize_Redundant(DM dm,PetscInt rank,PetscInt N) 372 { 373 DM_Redundant *red = (DM_Redundant*)dm->data; 374 PetscErrorCode ierr; 375 PetscMPIInt myrank; 376 377 PetscFunctionBegin; 378 ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)dm),&myrank);CHKERRQ(ierr); 379 red->rank = rank; 380 red->N = N; 381 red->n = (myrank == rank) ? N : 0; 382 PetscFunctionReturn(0); 383 } 384 385 #undef __FUNCT__ 386 #define __FUNCT__ "DMRedundantGetSize_Redundant" 387 static PetscErrorCode DMRedundantGetSize_Redundant(DM dm,PetscInt *rank,PetscInt *N) 388 { 389 DM_Redundant *red = (DM_Redundant*)dm->data; 390 391 PetscFunctionBegin; 392 if (rank) *rank = red->rank; 393 if (N) *N = red->N; 394 PetscFunctionReturn(0); 395 } 396 397 /*MC 398 DMREDUNDANT = "redundant" - A DM object that is used to manage data for a small set of dense globally coupled variables. 399 In the global representation of the vector the variables are all stored on a single MPI process (all the other MPI processes have 400 no variables) in the local representation all the variables are stored on ALL the MPI processes (because they are all needed for each 401 processes local computations). 402 403 This DM is generally used inside a DMCOMPOSITE object. For example, it may be used to store continuation parameters for a bifurcation problem. 404 405 406 Level: intermediate 407 408 .seealso: DMType, DMCOMPOSITE, DMCreateRedundant(), DMCreate(), DMRedundantSetSize(), DMRedundantGetSize() 409 M*/ 410 411 #undef __FUNCT__ 412 #define __FUNCT__ "DMCreate_Redundant" 413 PETSC_EXTERN PetscErrorCode DMCreate_Redundant(DM dm) 414 { 415 PetscErrorCode ierr; 416 DM_Redundant *red; 417 418 PetscFunctionBegin; 419 ierr = PetscNewLog(dm,&red);CHKERRQ(ierr); 420 dm->data = red; 421 422 ierr = PetscObjectChangeTypeName((PetscObject)dm,DMREDUNDANT);CHKERRQ(ierr); 423 424 dm->ops->setup = DMSetUp_Redundant; 425 dm->ops->view = DMView_Redundant; 426 dm->ops->createglobalvector = DMCreateGlobalVector_Redundant; 427 dm->ops->createlocalvector = DMCreateLocalVector_Redundant; 428 dm->ops->creatematrix = DMCreateMatrix_Redundant; 429 dm->ops->destroy = DMDestroy_Redundant; 430 dm->ops->globaltolocalbegin = DMGlobalToLocalBegin_Redundant; 431 dm->ops->globaltolocalend = DMGlobalToLocalEnd_Redundant; 432 dm->ops->localtoglobalbegin = DMLocalToGlobalBegin_Redundant; 433 dm->ops->localtoglobalend = DMLocalToGlobalEnd_Redundant; 434 dm->ops->refine = DMRefine_Redundant; 435 dm->ops->coarsen = DMCoarsen_Redundant; 436 dm->ops->createinterpolation= DMCreateInterpolation_Redundant; 437 dm->ops->getcoloring = DMCreateColoring_Redundant; 438 439 ierr = PetscObjectComposeFunction((PetscObject)dm,"DMRedundantSetSize_C",DMRedundantSetSize_Redundant);CHKERRQ(ierr); 440 ierr = PetscObjectComposeFunction((PetscObject)dm,"DMRedundantGetSize_C",DMRedundantGetSize_Redundant);CHKERRQ(ierr); 441 PetscFunctionReturn(0); 442 } 443 444 #undef __FUNCT__ 445 #define __FUNCT__ "DMRedundantCreate" 446 /*@C 447 DMRedundantCreate - Creates a DM object, used to manage data for dense globally coupled variables 448 449 Collective on MPI_Comm 450 451 Input Parameter: 452 + comm - the processors that will share the global vector 453 . rank - rank to own the redundant values 454 - N - total number of degrees of freedom 455 456 Output Parameters: 457 . red - the redundant DM 458 459 Level: advanced 460 461 .seealso DMDestroy(), DMCreateGlobalVector(), DMCreateMatrix(), DMCompositeAddDM(), DMREDUNDANT, DMSetType(), DMRedundantSetSize(), DMRedundantGetSize() 462 463 @*/ 464 PetscErrorCode DMRedundantCreate(MPI_Comm comm,PetscInt rank,PetscInt N,DM *dm) 465 { 466 PetscErrorCode ierr; 467 468 PetscFunctionBegin; 469 PetscValidPointer(dm,2); 470 ierr = DMCreate(comm,dm);CHKERRQ(ierr); 471 ierr = DMSetType(*dm,DMREDUNDANT);CHKERRQ(ierr); 472 ierr = DMRedundantSetSize(*dm,rank,N);CHKERRQ(ierr); 473 ierr = DMSetUp(*dm);CHKERRQ(ierr); 474 PetscFunctionReturn(0); 475 } 476