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