1*a7e14dcfSSatish Balay #include "tao.h" /*I "tao.h" I*/ 2*a7e14dcfSSatish Balay #include "petsc-private/matimpl.h" 3*a7e14dcfSSatish Balay #include "src/matrix/submatfree.h" 4*a7e14dcfSSatish Balay #include "tao_util.h" /*I "tao_util.h" I*/ 5*a7e14dcfSSatish Balay 6*a7e14dcfSSatish Balay 7*a7e14dcfSSatish Balay #undef __FUNCT__ 8*a7e14dcfSSatish Balay #define __FUNCT__ "VecWhichEqual" 9*a7e14dcfSSatish Balay /*@ 10*a7e14dcfSSatish Balay VecWhichEqual - Creates an index set containing the indices 11*a7e14dcfSSatish Balay where the vectors Vec1 and Vec2 have identical elements. 12*a7e14dcfSSatish Balay 13*a7e14dcfSSatish Balay Collective on S 14*a7e14dcfSSatish Balay 15*a7e14dcfSSatish Balay Input Parameters: 16*a7e14dcfSSatish Balay . Vec1, Vec2 - the two vectors to compare 17*a7e14dcfSSatish Balay 18*a7e14dcfSSatish Balay OutputParameter: 19*a7e14dcfSSatish Balay . S - The index set containing the indices i where vec1[i] == vec2[i] 20*a7e14dcfSSatish Balay 21*a7e14dcfSSatish Balay Level: advanced 22*a7e14dcfSSatish Balay @*/ 23*a7e14dcfSSatish Balay PetscErrorCode VecWhichEqual(Vec Vec1, Vec Vec2, IS * S) 24*a7e14dcfSSatish Balay { 25*a7e14dcfSSatish Balay PetscErrorCode ierr; 26*a7e14dcfSSatish Balay PetscInt i,n_same=0; 27*a7e14dcfSSatish Balay PetscInt n,low,high,low2,high2; 28*a7e14dcfSSatish Balay PetscInt *same; 29*a7e14dcfSSatish Balay PetscReal *v1,*v2; 30*a7e14dcfSSatish Balay MPI_Comm comm; 31*a7e14dcfSSatish Balay 32*a7e14dcfSSatish Balay PetscFunctionBegin; 33*a7e14dcfSSatish Balay PetscValidHeaderSpecific(Vec1,VEC_CLASSID,1); 34*a7e14dcfSSatish Balay PetscValidHeaderSpecific(Vec2,VEC_CLASSID,2); 35*a7e14dcfSSatish Balay PetscCheckSameComm(Vec1,1,Vec2,2); 36*a7e14dcfSSatish Balay 37*a7e14dcfSSatish Balay 38*a7e14dcfSSatish Balay ierr = VecGetOwnershipRange(Vec1, &low, &high); CHKERRQ(ierr); 39*a7e14dcfSSatish Balay ierr = VecGetOwnershipRange(Vec2, &low2, &high2); CHKERRQ(ierr); 40*a7e14dcfSSatish Balay 41*a7e14dcfSSatish Balay if ( low != low2 || high != high2 ) 42*a7e14dcfSSatish Balay SETERRQ(PETSC_COMM_SELF,1,"Vectors must be identically loaded over processors"); 43*a7e14dcfSSatish Balay 44*a7e14dcfSSatish Balay ierr = VecGetLocalSize(Vec1,&n); CHKERRQ(ierr); 45*a7e14dcfSSatish Balay 46*a7e14dcfSSatish Balay if (n>0){ 47*a7e14dcfSSatish Balay 48*a7e14dcfSSatish Balay if (Vec1 == Vec2){ 49*a7e14dcfSSatish Balay ierr = VecGetArray(Vec1,&v1); CHKERRQ(ierr); 50*a7e14dcfSSatish Balay v2=v1; 51*a7e14dcfSSatish Balay } else { 52*a7e14dcfSSatish Balay ierr = VecGetArray(Vec1,&v1); CHKERRQ(ierr); 53*a7e14dcfSSatish Balay ierr = VecGetArray(Vec2,&v2); CHKERRQ(ierr); 54*a7e14dcfSSatish Balay } 55*a7e14dcfSSatish Balay 56*a7e14dcfSSatish Balay ierr = PetscMalloc( n*sizeof(PetscInt),&same ); CHKERRQ(ierr); 57*a7e14dcfSSatish Balay 58*a7e14dcfSSatish Balay for (i=0; i<n; i++){ 59*a7e14dcfSSatish Balay if (v1[i] == v2[i]) {same[n_same]=low+i; n_same++;} 60*a7e14dcfSSatish Balay } 61*a7e14dcfSSatish Balay 62*a7e14dcfSSatish Balay if (Vec1 == Vec2){ 63*a7e14dcfSSatish Balay ierr = VecRestoreArray(Vec1,&v1); CHKERRQ(ierr); 64*a7e14dcfSSatish Balay } else { 65*a7e14dcfSSatish Balay ierr = VecRestoreArray(Vec1,&v1); CHKERRQ(ierr); 66*a7e14dcfSSatish Balay ierr = VecRestoreArray(Vec2,&v2); CHKERRQ(ierr); 67*a7e14dcfSSatish Balay } 68*a7e14dcfSSatish Balay 69*a7e14dcfSSatish Balay } else { 70*a7e14dcfSSatish Balay 71*a7e14dcfSSatish Balay n_same = 0; same=NULL; 72*a7e14dcfSSatish Balay 73*a7e14dcfSSatish Balay } 74*a7e14dcfSSatish Balay 75*a7e14dcfSSatish Balay ierr = PetscObjectGetComm((PetscObject)Vec1,&comm);CHKERRQ(ierr); 76*a7e14dcfSSatish Balay ierr = ISCreateGeneral(comm,n_same,same,PETSC_COPY_VALUES,S);CHKERRQ(ierr); 77*a7e14dcfSSatish Balay 78*a7e14dcfSSatish Balay if (same) { 79*a7e14dcfSSatish Balay ierr = PetscFree(same); CHKERRQ(ierr); 80*a7e14dcfSSatish Balay } 81*a7e14dcfSSatish Balay 82*a7e14dcfSSatish Balay PetscFunctionReturn(0); 83*a7e14dcfSSatish Balay } 84*a7e14dcfSSatish Balay 85*a7e14dcfSSatish Balay #undef __FUNCT__ 86*a7e14dcfSSatish Balay #define __FUNCT__ "VecWhichLessThan" 87*a7e14dcfSSatish Balay /*@ 88*a7e14dcfSSatish Balay VecWhichLessThan - Creates an index set containing the indices 89*a7e14dcfSSatish Balay where the vectors Vec1 < Vec2 90*a7e14dcfSSatish Balay 91*a7e14dcfSSatish Balay Collective on S 92*a7e14dcfSSatish Balay 93*a7e14dcfSSatish Balay Input Parameters: 94*a7e14dcfSSatish Balay . Vec1, Vec2 - the two vectors to compare 95*a7e14dcfSSatish Balay 96*a7e14dcfSSatish Balay OutputParameter: 97*a7e14dcfSSatish Balay . S - The index set containing the indices i where vec1[i] < vec2[i] 98*a7e14dcfSSatish Balay 99*a7e14dcfSSatish Balay Level: advanced 100*a7e14dcfSSatish Balay @*/ 101*a7e14dcfSSatish Balay PetscErrorCode VecWhichLessThan(Vec Vec1, Vec Vec2, IS * S) 102*a7e14dcfSSatish Balay { 103*a7e14dcfSSatish Balay int ierr; 104*a7e14dcfSSatish Balay PetscInt i; 105*a7e14dcfSSatish Balay PetscInt n,low,high,low2,high2,n_lt=0; 106*a7e14dcfSSatish Balay PetscInt *lt; 107*a7e14dcfSSatish Balay PetscReal *v1,*v2; 108*a7e14dcfSSatish Balay MPI_Comm comm; 109*a7e14dcfSSatish Balay 110*a7e14dcfSSatish Balay PetscFunctionBegin; 111*a7e14dcfSSatish Balay PetscValidHeaderSpecific(Vec1,VEC_CLASSID,1); 112*a7e14dcfSSatish Balay PetscValidHeaderSpecific(Vec2,VEC_CLASSID,2); 113*a7e14dcfSSatish Balay PetscCheckSameComm(Vec1,1,Vec2,2); 114*a7e14dcfSSatish Balay 115*a7e14dcfSSatish Balay ierr = VecGetOwnershipRange(Vec1, &low, &high); CHKERRQ(ierr); 116*a7e14dcfSSatish Balay ierr = VecGetOwnershipRange(Vec2, &low2, &high2); CHKERRQ(ierr); 117*a7e14dcfSSatish Balay 118*a7e14dcfSSatish Balay if ( low != low2 || high != high2 ) 119*a7e14dcfSSatish Balay SETERRQ(PETSC_COMM_SELF,1,"Vectors must be identically loaded over processors"); 120*a7e14dcfSSatish Balay 121*a7e14dcfSSatish Balay ierr = VecGetLocalSize(Vec1,&n); CHKERRQ(ierr); 122*a7e14dcfSSatish Balay 123*a7e14dcfSSatish Balay if (n>0){ 124*a7e14dcfSSatish Balay 125*a7e14dcfSSatish Balay if (Vec1 == Vec2){ 126*a7e14dcfSSatish Balay ierr = VecGetArray(Vec1,&v1); CHKERRQ(ierr); 127*a7e14dcfSSatish Balay v2=v1; 128*a7e14dcfSSatish Balay } else { 129*a7e14dcfSSatish Balay ierr = VecGetArray(Vec1,&v1); CHKERRQ(ierr); 130*a7e14dcfSSatish Balay ierr = VecGetArray(Vec2,&v2); CHKERRQ(ierr); 131*a7e14dcfSSatish Balay } 132*a7e14dcfSSatish Balay ierr = PetscMalloc( n*sizeof(PetscInt),< ); CHKERRQ(ierr); 133*a7e14dcfSSatish Balay 134*a7e14dcfSSatish Balay for (i=0; i<n; i++){ 135*a7e14dcfSSatish Balay if (v1[i] < v2[i]) {lt[n_lt]=low+i; n_lt++;} 136*a7e14dcfSSatish Balay } 137*a7e14dcfSSatish Balay 138*a7e14dcfSSatish Balay if (Vec1 == Vec2){ 139*a7e14dcfSSatish Balay ierr = VecRestoreArray(Vec1,&v1); CHKERRQ(ierr); 140*a7e14dcfSSatish Balay } else { 141*a7e14dcfSSatish Balay ierr = VecRestoreArray(Vec1,&v1); CHKERRQ(ierr); 142*a7e14dcfSSatish Balay ierr = VecRestoreArray(Vec2,&v2); CHKERRQ(ierr); 143*a7e14dcfSSatish Balay } 144*a7e14dcfSSatish Balay 145*a7e14dcfSSatish Balay } else { 146*a7e14dcfSSatish Balay n_lt=0; lt=NULL; 147*a7e14dcfSSatish Balay } 148*a7e14dcfSSatish Balay 149*a7e14dcfSSatish Balay ierr = PetscObjectGetComm((PetscObject)Vec1,&comm);CHKERRQ(ierr); 150*a7e14dcfSSatish Balay ierr = ISCreateGeneral(comm,n_lt,lt,PETSC_COPY_VALUES,S);CHKERRQ(ierr); 151*a7e14dcfSSatish Balay 152*a7e14dcfSSatish Balay if (lt) { 153*a7e14dcfSSatish Balay ierr = PetscFree(lt); CHKERRQ(ierr); 154*a7e14dcfSSatish Balay } 155*a7e14dcfSSatish Balay 156*a7e14dcfSSatish Balay PetscFunctionReturn(0); 157*a7e14dcfSSatish Balay } 158*a7e14dcfSSatish Balay 159*a7e14dcfSSatish Balay #undef __FUNCT__ 160*a7e14dcfSSatish Balay #define __FUNCT__ "VecWhichGreaterThan" 161*a7e14dcfSSatish Balay /*@ 162*a7e14dcfSSatish Balay VecWhichGreaterThan - Creates an index set containing the indices 163*a7e14dcfSSatish Balay where the vectors Vec1 > Vec2 164*a7e14dcfSSatish Balay 165*a7e14dcfSSatish Balay Collective on S 166*a7e14dcfSSatish Balay 167*a7e14dcfSSatish Balay Input Parameters: 168*a7e14dcfSSatish Balay . Vec1, Vec2 - the two vectors to compare 169*a7e14dcfSSatish Balay 170*a7e14dcfSSatish Balay OutputParameter: 171*a7e14dcfSSatish Balay . S - The index set containing the indices i where vec1[i] > vec2[i] 172*a7e14dcfSSatish Balay 173*a7e14dcfSSatish Balay Level: advanced 174*a7e14dcfSSatish Balay @*/ 175*a7e14dcfSSatish Balay PetscErrorCode VecWhichGreaterThan(Vec Vec1, Vec Vec2, IS * S) 176*a7e14dcfSSatish Balay { 177*a7e14dcfSSatish Balay int ierr; 178*a7e14dcfSSatish Balay PetscInt n,low,high,low2,high2,n_gt=0,i; 179*a7e14dcfSSatish Balay PetscInt *gt=NULL; 180*a7e14dcfSSatish Balay PetscReal *v1,*v2; 181*a7e14dcfSSatish Balay MPI_Comm comm; 182*a7e14dcfSSatish Balay 183*a7e14dcfSSatish Balay PetscFunctionBegin; 184*a7e14dcfSSatish Balay PetscValidHeaderSpecific(Vec1,VEC_CLASSID,1); 185*a7e14dcfSSatish Balay PetscValidHeaderSpecific(Vec2,VEC_CLASSID,2); 186*a7e14dcfSSatish Balay PetscCheckSameComm(Vec1,1,Vec2,2); 187*a7e14dcfSSatish Balay 188*a7e14dcfSSatish Balay ierr = VecGetOwnershipRange(Vec1, &low, &high); CHKERRQ(ierr); 189*a7e14dcfSSatish Balay ierr = VecGetOwnershipRange(Vec2, &low2, &high2); CHKERRQ(ierr); 190*a7e14dcfSSatish Balay 191*a7e14dcfSSatish Balay if ( low != low2 || high != high2 ) 192*a7e14dcfSSatish Balay SETERRQ(PETSC_COMM_SELF,1,"Vectors must be identically loaded over processors"); 193*a7e14dcfSSatish Balay 194*a7e14dcfSSatish Balay ierr = VecGetLocalSize(Vec1,&n); CHKERRQ(ierr); 195*a7e14dcfSSatish Balay 196*a7e14dcfSSatish Balay if (n>0){ 197*a7e14dcfSSatish Balay 198*a7e14dcfSSatish Balay if (Vec1 == Vec2){ 199*a7e14dcfSSatish Balay ierr = VecGetArray(Vec1,&v1); CHKERRQ(ierr); 200*a7e14dcfSSatish Balay v2=v1; 201*a7e14dcfSSatish Balay } else { 202*a7e14dcfSSatish Balay ierr = VecGetArray(Vec1,&v1); CHKERRQ(ierr); 203*a7e14dcfSSatish Balay ierr = VecGetArray(Vec2,&v2); CHKERRQ(ierr); 204*a7e14dcfSSatish Balay } 205*a7e14dcfSSatish Balay 206*a7e14dcfSSatish Balay ierr = PetscMalloc( n*sizeof(PetscInt), > ); CHKERRQ(ierr); 207*a7e14dcfSSatish Balay 208*a7e14dcfSSatish Balay for (i=0; i<n; i++){ 209*a7e14dcfSSatish Balay if (v1[i] > v2[i]) {gt[n_gt]=low+i; n_gt++;} 210*a7e14dcfSSatish Balay } 211*a7e14dcfSSatish Balay 212*a7e14dcfSSatish Balay if (Vec1 == Vec2){ 213*a7e14dcfSSatish Balay ierr = VecRestoreArray(Vec1,&v1); CHKERRQ(ierr); 214*a7e14dcfSSatish Balay } else { 215*a7e14dcfSSatish Balay ierr = VecRestoreArray(Vec1,&v1); CHKERRQ(ierr); 216*a7e14dcfSSatish Balay ierr = VecRestoreArray(Vec2,&v2); CHKERRQ(ierr); 217*a7e14dcfSSatish Balay } 218*a7e14dcfSSatish Balay 219*a7e14dcfSSatish Balay } else{ 220*a7e14dcfSSatish Balay 221*a7e14dcfSSatish Balay n_gt=0; gt=NULL; 222*a7e14dcfSSatish Balay 223*a7e14dcfSSatish Balay } 224*a7e14dcfSSatish Balay 225*a7e14dcfSSatish Balay ierr = PetscObjectGetComm((PetscObject)Vec1,&comm);CHKERRQ(ierr); 226*a7e14dcfSSatish Balay ierr = ISCreateGeneral(comm,n_gt,gt,PETSC_COPY_VALUES,S);CHKERRQ(ierr); 227*a7e14dcfSSatish Balay 228*a7e14dcfSSatish Balay if (gt) { 229*a7e14dcfSSatish Balay ierr = PetscFree(gt); CHKERRQ(ierr); 230*a7e14dcfSSatish Balay } 231*a7e14dcfSSatish Balay PetscFunctionReturn(0); 232*a7e14dcfSSatish Balay } 233*a7e14dcfSSatish Balay 234*a7e14dcfSSatish Balay #undef __FUNCT__ 235*a7e14dcfSSatish Balay #define __FUNCT__ "VecWhichBetween" 236*a7e14dcfSSatish Balay /*@ 237*a7e14dcfSSatish Balay VecWhichBetween - Creates an index set containing the indices 238*a7e14dcfSSatish Balay where VecLow < V < VecHigh 239*a7e14dcfSSatish Balay 240*a7e14dcfSSatish Balay Collective on S 241*a7e14dcfSSatish Balay 242*a7e14dcfSSatish Balay Input Parameters: 243*a7e14dcfSSatish Balay + VecLow - lower bound 244*a7e14dcfSSatish Balay . V - Vector to compare 245*a7e14dcfSSatish Balay - VecHigh - higher bound 246*a7e14dcfSSatish Balay 247*a7e14dcfSSatish Balay OutputParameter: 248*a7e14dcfSSatish Balay . S - The index set containing the indices i where veclow[i] < v[i] < vechigh[i] 249*a7e14dcfSSatish Balay 250*a7e14dcfSSatish Balay Level: advanced 251*a7e14dcfSSatish Balay @*/ 252*a7e14dcfSSatish Balay PetscErrorCode VecWhichBetween(Vec VecLow, Vec V, Vec VecHigh, IS *S) 253*a7e14dcfSSatish Balay { 254*a7e14dcfSSatish Balay 255*a7e14dcfSSatish Balay PetscErrorCode ierr; 256*a7e14dcfSSatish Balay PetscInt n,low,high,low2,high2,low3,high3,n_vm=0; 257*a7e14dcfSSatish Balay PetscInt *vm,i; 258*a7e14dcfSSatish Balay PetscReal *v1,*v2,*vmiddle; 259*a7e14dcfSSatish Balay MPI_Comm comm; 260*a7e14dcfSSatish Balay 261*a7e14dcfSSatish Balay PetscValidHeaderSpecific(V,VEC_CLASSID,2); 262*a7e14dcfSSatish Balay PetscCheckSameComm(V,2,VecLow,1); PetscCheckSameComm(V,2,VecHigh,3); 263*a7e14dcfSSatish Balay 264*a7e14dcfSSatish Balay ierr = VecGetOwnershipRange(VecLow, &low, &high); CHKERRQ(ierr); 265*a7e14dcfSSatish Balay ierr = VecGetOwnershipRange(VecHigh, &low2, &high2); CHKERRQ(ierr); 266*a7e14dcfSSatish Balay ierr = VecGetOwnershipRange(V, &low3, &high3); CHKERRQ(ierr); 267*a7e14dcfSSatish Balay 268*a7e14dcfSSatish Balay if ( low!=low2 || high!=high2 || low!=low3 || high!=high3 ) 269*a7e14dcfSSatish Balay SETERRQ(PETSC_COMM_SELF,1,"Vectors must be identically loaded over processors"); 270*a7e14dcfSSatish Balay 271*a7e14dcfSSatish Balay ierr = VecGetLocalSize(VecLow,&n); CHKERRQ(ierr); 272*a7e14dcfSSatish Balay 273*a7e14dcfSSatish Balay if (n>0){ 274*a7e14dcfSSatish Balay 275*a7e14dcfSSatish Balay ierr = VecGetArray(VecLow,&v1); CHKERRQ(ierr); 276*a7e14dcfSSatish Balay if (VecLow != VecHigh){ 277*a7e14dcfSSatish Balay ierr = VecGetArray(VecHigh,&v2); CHKERRQ(ierr); 278*a7e14dcfSSatish Balay } else { 279*a7e14dcfSSatish Balay v2=v1; 280*a7e14dcfSSatish Balay } 281*a7e14dcfSSatish Balay if ( V != VecLow && V != VecHigh){ 282*a7e14dcfSSatish Balay ierr = VecGetArray(V,&vmiddle); CHKERRQ(ierr); 283*a7e14dcfSSatish Balay } else if ( V==VecLow ){ 284*a7e14dcfSSatish Balay vmiddle=v1; 285*a7e14dcfSSatish Balay } else { 286*a7e14dcfSSatish Balay vmiddle =v2; 287*a7e14dcfSSatish Balay } 288*a7e14dcfSSatish Balay 289*a7e14dcfSSatish Balay ierr = PetscMalloc( n*sizeof(PetscInt), &vm ); CHKERRQ(ierr); 290*a7e14dcfSSatish Balay 291*a7e14dcfSSatish Balay for (i=0; i<n; i++){ 292*a7e14dcfSSatish Balay if (v1[i] < vmiddle[i] && vmiddle[i] < v2[i]) {vm[n_vm]=low+i; n_vm++;} 293*a7e14dcfSSatish Balay } 294*a7e14dcfSSatish Balay 295*a7e14dcfSSatish Balay ierr = VecRestoreArray(VecLow,&v1); CHKERRQ(ierr); 296*a7e14dcfSSatish Balay if (VecLow != VecHigh){ 297*a7e14dcfSSatish Balay ierr = VecRestoreArray(VecHigh,&v2); CHKERRQ(ierr); 298*a7e14dcfSSatish Balay } 299*a7e14dcfSSatish Balay if ( V != VecLow && V != VecHigh){ 300*a7e14dcfSSatish Balay ierr = VecRestoreArray(V,&vmiddle); CHKERRQ(ierr); 301*a7e14dcfSSatish Balay } 302*a7e14dcfSSatish Balay 303*a7e14dcfSSatish Balay } else { 304*a7e14dcfSSatish Balay 305*a7e14dcfSSatish Balay n_vm=0; vm=NULL; 306*a7e14dcfSSatish Balay 307*a7e14dcfSSatish Balay } 308*a7e14dcfSSatish Balay 309*a7e14dcfSSatish Balay ierr = PetscObjectGetComm((PetscObject)V,&comm);CHKERRQ(ierr); 310*a7e14dcfSSatish Balay ierr = ISCreateGeneral(comm,n_vm,vm,PETSC_COPY_VALUES,S);CHKERRQ(ierr); 311*a7e14dcfSSatish Balay 312*a7e14dcfSSatish Balay if (vm) { 313*a7e14dcfSSatish Balay ierr = PetscFree(vm); CHKERRQ(ierr); 314*a7e14dcfSSatish Balay } 315*a7e14dcfSSatish Balay 316*a7e14dcfSSatish Balay PetscFunctionReturn(0); 317*a7e14dcfSSatish Balay } 318*a7e14dcfSSatish Balay 319*a7e14dcfSSatish Balay 320*a7e14dcfSSatish Balay #undef __FUNCT__ 321*a7e14dcfSSatish Balay #define __FUNCT__ "VecWhichBetweenOrEqual" 322*a7e14dcfSSatish Balay /*@ 323*a7e14dcfSSatish Balay VecWhichBetweenOrEqual - Creates an index set containing the indices 324*a7e14dcfSSatish Balay where VecLow <= V <= VecHigh 325*a7e14dcfSSatish Balay 326*a7e14dcfSSatish Balay Collective on S 327*a7e14dcfSSatish Balay 328*a7e14dcfSSatish Balay Input Parameters: 329*a7e14dcfSSatish Balay + VecLow - lower bound 330*a7e14dcfSSatish Balay . V - Vector to compare 331*a7e14dcfSSatish Balay - VecHigh - higher bound 332*a7e14dcfSSatish Balay 333*a7e14dcfSSatish Balay OutputParameter: 334*a7e14dcfSSatish Balay . S - The index set containing the indices i where veclow[i] <= v[i] <= vechigh[i] 335*a7e14dcfSSatish Balay 336*a7e14dcfSSatish Balay Level: advanced 337*a7e14dcfSSatish Balay @*/ 338*a7e14dcfSSatish Balay 339*a7e14dcfSSatish Balay PetscErrorCode VecWhichBetweenOrEqual(Vec VecLow, Vec V, Vec VecHigh, IS * S) 340*a7e14dcfSSatish Balay { 341*a7e14dcfSSatish Balay PetscErrorCode ierr; 342*a7e14dcfSSatish Balay PetscInt n,low,high,low2,high2,low3,high3,n_vm=0,i; 343*a7e14dcfSSatish Balay PetscInt *vm; 344*a7e14dcfSSatish Balay PetscReal *v1,*v2,*vmiddle; 345*a7e14dcfSSatish Balay MPI_Comm comm; 346*a7e14dcfSSatish Balay 347*a7e14dcfSSatish Balay PetscValidHeaderSpecific(V,VEC_CLASSID,2); 348*a7e14dcfSSatish Balay PetscCheckSameComm(V,2,VecLow,1); PetscCheckSameComm(V,2,VecHigh,3); 349*a7e14dcfSSatish Balay 350*a7e14dcfSSatish Balay ierr = VecGetOwnershipRange(VecLow, &low, &high); CHKERRQ(ierr); 351*a7e14dcfSSatish Balay ierr = VecGetOwnershipRange(VecHigh, &low2, &high2); CHKERRQ(ierr); 352*a7e14dcfSSatish Balay ierr = VecGetOwnershipRange(V, &low3, &high3); CHKERRQ(ierr); 353*a7e14dcfSSatish Balay 354*a7e14dcfSSatish Balay if ( low!=low2 || high!=high2 || low!=low3 || high!=high3 ) 355*a7e14dcfSSatish Balay SETERRQ(PETSC_COMM_SELF,1,"Vectors must be identically loaded over processors"); 356*a7e14dcfSSatish Balay 357*a7e14dcfSSatish Balay ierr = VecGetLocalSize(VecLow,&n); CHKERRQ(ierr); 358*a7e14dcfSSatish Balay 359*a7e14dcfSSatish Balay if (n>0){ 360*a7e14dcfSSatish Balay 361*a7e14dcfSSatish Balay ierr = VecGetArray(VecLow,&v1); CHKERRQ(ierr); 362*a7e14dcfSSatish Balay if (VecLow != VecHigh){ 363*a7e14dcfSSatish Balay ierr = VecGetArray(VecHigh,&v2); CHKERRQ(ierr); 364*a7e14dcfSSatish Balay } else { 365*a7e14dcfSSatish Balay v2=v1; 366*a7e14dcfSSatish Balay } 367*a7e14dcfSSatish Balay if ( V != VecLow && V != VecHigh){ 368*a7e14dcfSSatish Balay ierr = VecGetArray(V,&vmiddle); CHKERRQ(ierr); 369*a7e14dcfSSatish Balay } else if ( V==VecLow ){ 370*a7e14dcfSSatish Balay vmiddle=v1; 371*a7e14dcfSSatish Balay } else { 372*a7e14dcfSSatish Balay vmiddle =v2; 373*a7e14dcfSSatish Balay } 374*a7e14dcfSSatish Balay 375*a7e14dcfSSatish Balay ierr = PetscMalloc( n*sizeof(PetscInt), &vm ); CHKERRQ(ierr); 376*a7e14dcfSSatish Balay 377*a7e14dcfSSatish Balay for (i=0; i<n; i++){ 378*a7e14dcfSSatish Balay if (v1[i] <= vmiddle[i] && vmiddle[i] <= v2[i]) {vm[n_vm]=low+i; n_vm++;} 379*a7e14dcfSSatish Balay } 380*a7e14dcfSSatish Balay 381*a7e14dcfSSatish Balay ierr = VecRestoreArray(VecLow,&v1); CHKERRQ(ierr); 382*a7e14dcfSSatish Balay if (VecLow != VecHigh){ 383*a7e14dcfSSatish Balay ierr = VecRestoreArray(VecHigh,&v2); CHKERRQ(ierr); 384*a7e14dcfSSatish Balay } 385*a7e14dcfSSatish Balay if ( V != VecLow && V != VecHigh){ 386*a7e14dcfSSatish Balay ierr = VecRestoreArray(V,&vmiddle); CHKERRQ(ierr); 387*a7e14dcfSSatish Balay } 388*a7e14dcfSSatish Balay 389*a7e14dcfSSatish Balay } else { 390*a7e14dcfSSatish Balay 391*a7e14dcfSSatish Balay n_vm=0; vm=NULL; 392*a7e14dcfSSatish Balay 393*a7e14dcfSSatish Balay } 394*a7e14dcfSSatish Balay 395*a7e14dcfSSatish Balay ierr = PetscObjectGetComm((PetscObject)V,&comm);CHKERRQ(ierr); 396*a7e14dcfSSatish Balay ierr = ISCreateGeneral(comm,n_vm,vm,PETSC_COPY_VALUES,S);CHKERRQ(ierr); 397*a7e14dcfSSatish Balay 398*a7e14dcfSSatish Balay if (vm) { 399*a7e14dcfSSatish Balay ierr = PetscFree(vm); CHKERRQ(ierr); 400*a7e14dcfSSatish Balay } 401*a7e14dcfSSatish Balay 402*a7e14dcfSSatish Balay PetscFunctionReturn(0); 403*a7e14dcfSSatish Balay } 404*a7e14dcfSSatish Balay 405*a7e14dcfSSatish Balay 406*a7e14dcfSSatish Balay #undef __FUNCT__ 407*a7e14dcfSSatish Balay #define __FUNCT__ "VecGetSubVec" 408*a7e14dcfSSatish Balay /*@ 409*a7e14dcfSSatish Balay VecGetSubVec - Gets a subvector using the IS 410*a7e14dcfSSatish Balay 411*a7e14dcfSSatish Balay Input Parameters: 412*a7e14dcfSSatish Balay + vfull - the full matrix 413*a7e14dcfSSatish Balay . is - the index set for the subvector 414*a7e14dcfSSatish Balay . reduced_type - the method TAO is using for subsetting (TAO_SUBSET_SUBVEC, TAO_SUBSET_MASK, TAO_SUBSET_MATRIXFREE) 415*a7e14dcfSSatish Balay - maskvalue - the value to set the unused vector elements to (for TAO_SUBSET_MASK or TAO_SUBSET_MATRIXFREE) 416*a7e14dcfSSatish Balay 417*a7e14dcfSSatish Balay 418*a7e14dcfSSatish Balay Output Parameters: 419*a7e14dcfSSatish Balay . vreduced - the subvector 420*a7e14dcfSSatish Balay 421*a7e14dcfSSatish Balay Note: 422*a7e14dcfSSatish Balay maskvalue should usually be 0.0, unless a pointwise divide will be used. 423*a7e14dcfSSatish Balay @*/ 424*a7e14dcfSSatish Balay PetscErrorCode VecGetSubVec(Vec vfull, IS is, PetscInt reduced_type, PetscReal maskvalue, Vec *vreduced) 425*a7e14dcfSSatish Balay { 426*a7e14dcfSSatish Balay PetscErrorCode ierr; 427*a7e14dcfSSatish Balay PetscInt nfull,nreduced,nreduced_local,rlow,rhigh,flow,fhigh; 428*a7e14dcfSSatish Balay PetscInt i,nlocal; 429*a7e14dcfSSatish Balay PetscReal *fv,*rv; 430*a7e14dcfSSatish Balay const PetscInt *s; 431*a7e14dcfSSatish Balay IS ident; 432*a7e14dcfSSatish Balay VecType vtype; 433*a7e14dcfSSatish Balay VecScatter scatter; 434*a7e14dcfSSatish Balay MPI_Comm comm; 435*a7e14dcfSSatish Balay 436*a7e14dcfSSatish Balay PetscFunctionBegin; 437*a7e14dcfSSatish Balay PetscValidHeaderSpecific(vfull,VEC_CLASSID,1); 438*a7e14dcfSSatish Balay PetscValidHeaderSpecific(is,IS_CLASSID,2); 439*a7e14dcfSSatish Balay 440*a7e14dcfSSatish Balay ierr = VecGetSize(vfull, &nfull); CHKERRQ(ierr); 441*a7e14dcfSSatish Balay ierr = ISGetSize(is, &nreduced); CHKERRQ(ierr); 442*a7e14dcfSSatish Balay 443*a7e14dcfSSatish Balay if (nreduced == nfull) { 444*a7e14dcfSSatish Balay 445*a7e14dcfSSatish Balay ierr = VecDestroy(vreduced); CHKERRQ(ierr); 446*a7e14dcfSSatish Balay ierr = VecDuplicate(vfull,vreduced); CHKERRQ(ierr); 447*a7e14dcfSSatish Balay ierr = VecCopy(vfull,*vreduced); CHKERRQ(ierr); 448*a7e14dcfSSatish Balay 449*a7e14dcfSSatish Balay } else { 450*a7e14dcfSSatish Balay 451*a7e14dcfSSatish Balay switch (reduced_type) { 452*a7e14dcfSSatish Balay case TAO_SUBSET_SUBVEC: 453*a7e14dcfSSatish Balay ierr = VecGetType(vfull,&vtype); CHKERRQ(ierr); 454*a7e14dcfSSatish Balay ierr = VecGetOwnershipRange(vfull,&flow,&fhigh); CHKERRQ(ierr); 455*a7e14dcfSSatish Balay ierr = ISGetLocalSize(is,&nreduced_local); CHKERRQ(ierr); 456*a7e14dcfSSatish Balay ierr = PetscObjectGetComm((PetscObject)vfull,&comm); CHKERRQ(ierr); 457*a7e14dcfSSatish Balay if (*vreduced) { 458*a7e14dcfSSatish Balay ierr = VecDestroy(vreduced); CHKERRQ(ierr); 459*a7e14dcfSSatish Balay } 460*a7e14dcfSSatish Balay ierr = VecCreate(comm,vreduced); CHKERRQ(ierr); 461*a7e14dcfSSatish Balay ierr = VecSetType(*vreduced,vtype); CHKERRQ(ierr); 462*a7e14dcfSSatish Balay 463*a7e14dcfSSatish Balay 464*a7e14dcfSSatish Balay ierr = VecSetSizes(*vreduced,nreduced_local,nreduced); CHKERRQ(ierr); 465*a7e14dcfSSatish Balay 466*a7e14dcfSSatish Balay ierr = VecGetOwnershipRange(*vreduced,&rlow,&rhigh); CHKERRQ(ierr); 467*a7e14dcfSSatish Balay 468*a7e14dcfSSatish Balay ierr = ISCreateStride(comm,nreduced_local,rlow,1,&ident); CHKERRQ(ierr); 469*a7e14dcfSSatish Balay ierr = VecScatterCreate(vfull,is,*vreduced,ident,&scatter); CHKERRQ(ierr); 470*a7e14dcfSSatish Balay ierr = VecScatterBegin(scatter,vfull,*vreduced,INSERT_VALUES,SCATTER_FORWARD); CHKERRQ(ierr); 471*a7e14dcfSSatish Balay ierr = VecScatterEnd(scatter,vfull,*vreduced,INSERT_VALUES,SCATTER_FORWARD); CHKERRQ(ierr); 472*a7e14dcfSSatish Balay ierr = VecScatterDestroy(&scatter); CHKERRQ(ierr); 473*a7e14dcfSSatish Balay ierr = ISDestroy(&ident); CHKERRQ(ierr); 474*a7e14dcfSSatish Balay break; 475*a7e14dcfSSatish Balay 476*a7e14dcfSSatish Balay case TAO_SUBSET_MASK: 477*a7e14dcfSSatish Balay case TAO_SUBSET_MATRIXFREE: 478*a7e14dcfSSatish Balay /* vr[i] = vf[i] if i in is 479*a7e14dcfSSatish Balay vr[i] = 0 otherwise */ 480*a7e14dcfSSatish Balay if (*vreduced == PETSC_NULL) { 481*a7e14dcfSSatish Balay ierr = VecDuplicate(vfull,vreduced); CHKERRQ(ierr); 482*a7e14dcfSSatish Balay } 483*a7e14dcfSSatish Balay 484*a7e14dcfSSatish Balay ierr = VecSet(*vreduced,maskvalue); CHKERRQ(ierr); 485*a7e14dcfSSatish Balay ierr = ISGetLocalSize(is,&nlocal); CHKERRQ(ierr); 486*a7e14dcfSSatish Balay ierr = VecGetOwnershipRange(vfull,&flow,&fhigh); CHKERRQ(ierr); 487*a7e14dcfSSatish Balay ierr = VecGetArray(vfull,&fv); CHKERRQ(ierr); 488*a7e14dcfSSatish Balay ierr = VecGetArray(*vreduced,&rv); CHKERRQ(ierr); 489*a7e14dcfSSatish Balay ierr = ISGetIndices(is,&s); CHKERRQ(ierr); 490*a7e14dcfSSatish Balay if (nlocal > (fhigh-flow)) { 491*a7e14dcfSSatish Balay SETERRQ2(PETSC_COMM_WORLD,1,"IS local size %d > Vec local size %d",nlocal,fhigh-flow); 492*a7e14dcfSSatish Balay } 493*a7e14dcfSSatish Balay for (i=0;i<nlocal;i++) { 494*a7e14dcfSSatish Balay rv[s[i]-flow] = fv[s[i]-flow]; 495*a7e14dcfSSatish Balay } 496*a7e14dcfSSatish Balay ierr = ISRestoreIndices(is,&s); CHKERRQ(ierr); 497*a7e14dcfSSatish Balay ierr = VecRestoreArray(vfull,&fv); CHKERRQ(ierr); 498*a7e14dcfSSatish Balay ierr = VecRestoreArray(*vreduced,&rv); CHKERRQ(ierr); 499*a7e14dcfSSatish Balay break; 500*a7e14dcfSSatish Balay } 501*a7e14dcfSSatish Balay } 502*a7e14dcfSSatish Balay PetscFunctionReturn(0); 503*a7e14dcfSSatish Balay 504*a7e14dcfSSatish Balay } 505*a7e14dcfSSatish Balay 506*a7e14dcfSSatish Balay #undef __FUNCT__ 507*a7e14dcfSSatish Balay #define __FUNCT__ "VecReducedXPY" 508*a7e14dcfSSatish Balay /*@ 509*a7e14dcfSSatish Balay VecReducedXPY - Adds a reduced vector to the appropriate elements of a full-space vector. 510*a7e14dcfSSatish Balay 511*a7e14dcfSSatish Balay Input Parameters: 512*a7e14dcfSSatish Balay + vfull - the full-space vector 513*a7e14dcfSSatish Balay . vreduced - the reduced-space vector 514*a7e14dcfSSatish Balay - is - the index set for the reduced space 515*a7e14dcfSSatish Balay 516*a7e14dcfSSatish Balay Output Parameters: 517*a7e14dcfSSatish Balay . vfull - the sum of the full-space vector and reduced-space vector 518*a7e14dcfSSatish Balay @*/ 519*a7e14dcfSSatish Balay PetscErrorCode VecReducedXPY(Vec vfull, Vec vreduced, IS is) 520*a7e14dcfSSatish Balay { 521*a7e14dcfSSatish Balay VecScatter scatter; 522*a7e14dcfSSatish Balay IS ident; 523*a7e14dcfSSatish Balay PetscInt nfull,nreduced,rlow,rhigh; 524*a7e14dcfSSatish Balay MPI_Comm comm; 525*a7e14dcfSSatish Balay PetscErrorCode ierr; 526*a7e14dcfSSatish Balay 527*a7e14dcfSSatish Balay PetscFunctionBegin; 528*a7e14dcfSSatish Balay PetscValidHeaderSpecific(vfull,VEC_CLASSID,1); 529*a7e14dcfSSatish Balay PetscValidHeaderSpecific(vreduced,VEC_CLASSID,2); 530*a7e14dcfSSatish Balay PetscValidHeaderSpecific(is,IS_CLASSID,3); 531*a7e14dcfSSatish Balay ierr = VecGetSize(vfull,&nfull); CHKERRQ(ierr); 532*a7e14dcfSSatish Balay ierr = VecGetSize(vreduced,&nreduced); CHKERRQ(ierr); 533*a7e14dcfSSatish Balay 534*a7e14dcfSSatish Balay if (nfull == nreduced) { /* Also takes care of masked vectors */ 535*a7e14dcfSSatish Balay ierr = VecAXPY(vfull,1.0,vreduced); CHKERRQ(ierr); 536*a7e14dcfSSatish Balay } else { 537*a7e14dcfSSatish Balay ierr = PetscObjectGetComm((PetscObject)vfull,&comm); CHKERRQ(ierr); 538*a7e14dcfSSatish Balay ierr = VecGetOwnershipRange(vreduced,&rlow,&rhigh); CHKERRQ(ierr); 539*a7e14dcfSSatish Balay ierr = ISCreateStride(comm,rhigh-rlow,rlow,1,&ident); CHKERRQ(ierr); 540*a7e14dcfSSatish Balay ierr = VecScatterCreate(vreduced,ident,vfull,is,&scatter); CHKERRQ(ierr); 541*a7e14dcfSSatish Balay ierr = VecScatterBegin(scatter,vreduced,vfull,ADD_VALUES,SCATTER_FORWARD); CHKERRQ(ierr); 542*a7e14dcfSSatish Balay ierr = VecScatterEnd(scatter,vreduced,vfull,ADD_VALUES,SCATTER_FORWARD); CHKERRQ(ierr); 543*a7e14dcfSSatish Balay ierr = VecScatterDestroy(&scatter); CHKERRQ(ierr); 544*a7e14dcfSSatish Balay ierr = ISDestroy(&ident); CHKERRQ(ierr); 545*a7e14dcfSSatish Balay } 546*a7e14dcfSSatish Balay 547*a7e14dcfSSatish Balay PetscFunctionReturn(0); 548*a7e14dcfSSatish Balay } 549*a7e14dcfSSatish Balay 550*a7e14dcfSSatish Balay 551*a7e14dcfSSatish Balay #undef __FUNCT__ 552*a7e14dcfSSatish Balay #define __FUNCT__ "ISCreateComplement" 553*a7e14dcfSSatish Balay /*@ 554*a7e14dcfSSatish Balay ISCreateComplement - Creates the complement of the the index set 555*a7e14dcfSSatish Balay 556*a7e14dcfSSatish Balay Collective on IS 557*a7e14dcfSSatish Balay 558*a7e14dcfSSatish Balay Input Parameter: 559*a7e14dcfSSatish Balay + S - a PETSc IS 560*a7e14dcfSSatish Balay - V - the reference vector space 561*a7e14dcfSSatish Balay 562*a7e14dcfSSatish Balay Output Parameter: 563*a7e14dcfSSatish Balay . T - the complement of S 564*a7e14dcfSSatish Balay 565*a7e14dcfSSatish Balay 566*a7e14dcfSSatish Balay .seealso ISCreateGeneral() 567*a7e14dcfSSatish Balay 568*a7e14dcfSSatish Balay Level: advanced 569*a7e14dcfSSatish Balay @*/ 570*a7e14dcfSSatish Balay PetscErrorCode ISCreateComplement(IS S, Vec V, IS *T){ 571*a7e14dcfSSatish Balay PetscErrorCode ierr; 572*a7e14dcfSSatish Balay PetscInt i,nis,nloc,high,low,n=0; 573*a7e14dcfSSatish Balay const PetscInt *s; 574*a7e14dcfSSatish Balay PetscInt *tt,*ss; 575*a7e14dcfSSatish Balay MPI_Comm comm; 576*a7e14dcfSSatish Balay 577*a7e14dcfSSatish Balay PetscFunctionBegin; 578*a7e14dcfSSatish Balay PetscValidHeaderSpecific(S,IS_CLASSID,1); 579*a7e14dcfSSatish Balay PetscValidHeaderSpecific(V,VEC_CLASSID,2); 580*a7e14dcfSSatish Balay 581*a7e14dcfSSatish Balay ierr = VecGetOwnershipRange(V,&low,&high); CHKERRQ(ierr); 582*a7e14dcfSSatish Balay ierr = VecGetLocalSize(V,&nloc); CHKERRQ(ierr); 583*a7e14dcfSSatish Balay ierr = ISGetLocalSize(S,&nis); CHKERRQ(ierr); 584*a7e14dcfSSatish Balay ierr = ISGetIndices(S, &s); CHKERRQ(ierr); 585*a7e14dcfSSatish Balay ierr = PetscMalloc( nloc*sizeof(PetscInt),&tt ); CHKERRQ(ierr); 586*a7e14dcfSSatish Balay ierr = PetscMalloc( nloc*sizeof(PetscInt),&ss ); CHKERRQ(ierr); 587*a7e14dcfSSatish Balay 588*a7e14dcfSSatish Balay for (i=low; i<high; i++){ tt[i-low]=i; } 589*a7e14dcfSSatish Balay 590*a7e14dcfSSatish Balay for (i=0; i<nis; i++){ tt[s[i]-low] = -2; } 591*a7e14dcfSSatish Balay 592*a7e14dcfSSatish Balay for (i=0; i<nloc; i++){ 593*a7e14dcfSSatish Balay if (tt[i]>-1){ ss[n]=tt[i]; n++; } 594*a7e14dcfSSatish Balay } 595*a7e14dcfSSatish Balay 596*a7e14dcfSSatish Balay ierr = ISRestoreIndices(S, &s); CHKERRQ(ierr); 597*a7e14dcfSSatish Balay 598*a7e14dcfSSatish Balay ierr = PetscObjectGetComm((PetscObject)S,&comm);CHKERRQ(ierr); 599*a7e14dcfSSatish Balay ierr = ISCreateGeneral(comm,n,ss,PETSC_COPY_VALUES,T);CHKERRQ(ierr); 600*a7e14dcfSSatish Balay 601*a7e14dcfSSatish Balay if (tt) { 602*a7e14dcfSSatish Balay ierr = PetscFree(tt); CHKERRQ(ierr); 603*a7e14dcfSSatish Balay } 604*a7e14dcfSSatish Balay if (ss) { 605*a7e14dcfSSatish Balay ierr = PetscFree(ss); CHKERRQ(ierr); 606*a7e14dcfSSatish Balay } 607*a7e14dcfSSatish Balay 608*a7e14dcfSSatish Balay PetscFunctionReturn(0); 609*a7e14dcfSSatish Balay } 610*a7e14dcfSSatish Balay 611*a7e14dcfSSatish Balay #undef __FUNCT__ 612*a7e14dcfSSatish Balay #define __FUNCT__ "VecISSetToConstant" 613*a7e14dcfSSatish Balay /*@ 614*a7e14dcfSSatish Balay VecISSetToConstant - Sets the elements of a vector, specified by an index set, to a constant 615*a7e14dcfSSatish Balay 616*a7e14dcfSSatish Balay Input Parameter: 617*a7e14dcfSSatish Balay + S - a PETSc IS 618*a7e14dcfSSatish Balay . c - the constant 619*a7e14dcfSSatish Balay - V - a Vec 620*a7e14dcfSSatish Balay 621*a7e14dcfSSatish Balay .seealso VecSet() 622*a7e14dcfSSatish Balay 623*a7e14dcfSSatish Balay Level: advanced 624*a7e14dcfSSatish Balay @*/ 625*a7e14dcfSSatish Balay PetscErrorCode VecISSetToConstant(IS S, PetscReal c, Vec V){ 626*a7e14dcfSSatish Balay PetscErrorCode ierr; 627*a7e14dcfSSatish Balay PetscInt nloc,low,high,i; 628*a7e14dcfSSatish Balay const PetscInt *s; 629*a7e14dcfSSatish Balay PetscReal *v; 630*a7e14dcfSSatish Balay PetscFunctionBegin; 631*a7e14dcfSSatish Balay PetscValidHeaderSpecific(V,VEC_CLASSID,3); 632*a7e14dcfSSatish Balay PetscValidHeaderSpecific(S,IS_CLASSID,1); 633*a7e14dcfSSatish Balay PetscValidType(V,3); 634*a7e14dcfSSatish Balay PetscCheckSameComm(V,3,S,1); 635*a7e14dcfSSatish Balay 636*a7e14dcfSSatish Balay ierr = VecGetOwnershipRange(V, &low, &high); CHKERRQ(ierr); 637*a7e14dcfSSatish Balay ierr = ISGetLocalSize(S,&nloc);CHKERRQ(ierr); 638*a7e14dcfSSatish Balay 639*a7e14dcfSSatish Balay ierr = ISGetIndices(S, &s); CHKERRQ(ierr); 640*a7e14dcfSSatish Balay ierr = VecGetArray(V,&v); CHKERRQ(ierr); 641*a7e14dcfSSatish Balay for (i=0; i<nloc; i++){ 642*a7e14dcfSSatish Balay v[s[i]-low] = c; 643*a7e14dcfSSatish Balay } 644*a7e14dcfSSatish Balay 645*a7e14dcfSSatish Balay ierr = ISRestoreIndices(S, &s); CHKERRQ(ierr); 646*a7e14dcfSSatish Balay ierr = VecRestoreArray(V,&v); CHKERRQ(ierr); 647*a7e14dcfSSatish Balay 648*a7e14dcfSSatish Balay PetscFunctionReturn(0); 649*a7e14dcfSSatish Balay 650*a7e14dcfSSatish Balay } 651*a7e14dcfSSatish Balay 652*a7e14dcfSSatish Balay #undef __FUNCT__ 653*a7e14dcfSSatish Balay #define __FUNCT__ "MatGetSubMat" 654*a7e14dcfSSatish Balay /*@ 655*a7e14dcfSSatish Balay MatGetSubMat - Gets a submatrix using the IS 656*a7e14dcfSSatish Balay 657*a7e14dcfSSatish Balay Input Parameters: 658*a7e14dcfSSatish Balay + M - the full matrix (n x n) 659*a7e14dcfSSatish Balay . is - the index set for the submatrix (both row and column index sets need to be the same) 660*a7e14dcfSSatish Balay . v1 - work vector of dimension n, needed for TAO_SUBSET_MASK option 661*a7e14dcfSSatish Balay - subset_type - the method TAO is using for subsetting (TAO_SUBSET_SUBVEC, TAO_SUBSET_MASK, 662*a7e14dcfSSatish Balay TAO_SUBSET_MATRIXFREE) 663*a7e14dcfSSatish Balay 664*a7e14dcfSSatish Balay Output Parameters: 665*a7e14dcfSSatish Balay . Msub - the submatrix 666*a7e14dcfSSatish Balay @*/ 667*a7e14dcfSSatish Balay PetscErrorCode MatGetSubMat(Mat M, IS is, Vec v1, TaoSubsetType subset_type, Mat *Msub) 668*a7e14dcfSSatish Balay { 669*a7e14dcfSSatish Balay PetscErrorCode ierr; 670*a7e14dcfSSatish Balay IS iscomp; 671*a7e14dcfSSatish Balay PetscBool flg; 672*a7e14dcfSSatish Balay PetscFunctionBegin; 673*a7e14dcfSSatish Balay PetscValidHeaderSpecific(M,MAT_CLASSID,1); 674*a7e14dcfSSatish Balay PetscValidHeaderSpecific(is,IS_CLASSID,2); 675*a7e14dcfSSatish Balay if (*Msub) { 676*a7e14dcfSSatish Balay ierr = MatDestroy(Msub); CHKERRQ(ierr); 677*a7e14dcfSSatish Balay } 678*a7e14dcfSSatish Balay switch (subset_type) { 679*a7e14dcfSSatish Balay case TAO_SUBSET_SUBVEC: 680*a7e14dcfSSatish Balay ierr = MatGetSubMatrix(M, is, is, MAT_INITIAL_MATRIX, Msub); CHKERRQ(ierr); 681*a7e14dcfSSatish Balay break; 682*a7e14dcfSSatish Balay 683*a7e14dcfSSatish Balay case TAO_SUBSET_MASK: 684*a7e14dcfSSatish Balay /* Get Reduced Hessian 685*a7e14dcfSSatish Balay Msub[i,j] = M[i,j] if i,j in Free_Local or i==j 686*a7e14dcfSSatish Balay Msub[i,j] = 0 if i!=j and i or j not in Free_Local 687*a7e14dcfSSatish Balay */ 688*a7e14dcfSSatish Balay ierr = PetscOptionsBool("-different_submatrix","use separate hessian matrix when computing submatrices","TaoSubsetType",PETSC_FALSE,&flg,PETSC_NULL); 689*a7e14dcfSSatish Balay if (flg == PETSC_TRUE) { 690*a7e14dcfSSatish Balay ierr = MatDuplicate(M, MAT_COPY_VALUES, Msub); CHKERRQ(ierr); 691*a7e14dcfSSatish Balay } else { 692*a7e14dcfSSatish Balay /* Act on hessian directly (default) */ 693*a7e14dcfSSatish Balay ierr = PetscObjectReference((PetscObject)M); CHKERRQ(ierr); 694*a7e14dcfSSatish Balay *Msub = M; 695*a7e14dcfSSatish Balay } 696*a7e14dcfSSatish Balay /* Save the diagonal to temporary vector */ 697*a7e14dcfSSatish Balay ierr = MatGetDiagonal(*Msub,v1); CHKERRQ(ierr); 698*a7e14dcfSSatish Balay 699*a7e14dcfSSatish Balay /* Zero out rows and columns */ 700*a7e14dcfSSatish Balay ierr = ISCreateComplement(is,v1,&iscomp); CHKERRQ(ierr); 701*a7e14dcfSSatish Balay 702*a7e14dcfSSatish Balay /* Use v1 instead of 0 here because of PETSc bug */ 703*a7e14dcfSSatish Balay ierr = MatZeroRowsColumnsIS(*Msub,iscomp,1.0,v1,v1); CHKERRQ(ierr); 704*a7e14dcfSSatish Balay 705*a7e14dcfSSatish Balay ierr = ISDestroy(&iscomp); CHKERRQ(ierr); 706*a7e14dcfSSatish Balay break; 707*a7e14dcfSSatish Balay case TAO_SUBSET_MATRIXFREE: 708*a7e14dcfSSatish Balay ierr = ISCreateComplement(is,v1,&iscomp); CHKERRQ(ierr); 709*a7e14dcfSSatish Balay ierr = MatCreateSubMatrixFree(M,iscomp,iscomp,Msub); CHKERRQ(ierr); 710*a7e14dcfSSatish Balay ierr = ISDestroy(&iscomp); CHKERRQ(ierr); 711*a7e14dcfSSatish Balay break; 712*a7e14dcfSSatish Balay } 713*a7e14dcfSSatish Balay PetscFunctionReturn(0); 714*a7e14dcfSSatish Balay } 715