Actual source code: isutil.c
petsc-dev 2014-02-02
1: #include <petsctao.h> /*I "petsctao.h" I*/
2: #include <petsc-private/matimpl.h>
3: #include <../src/tao/matrix/submatfree.h>
7: /*@
8: VecGetSubVec - Gets a subvector using the IS
10: Input Parameters:
11: + vfull - the full matrix
12: . is - the index set for the subvector
13: . reduced_type - the method TAO is using for subsetting (TAO_SUBSET_SUBVEC, TAO_SUBSET_MASK, TAO_SUBSET_MATRIXFREE)
14: - maskvalue - the value to set the unused vector elements to (for TAO_SUBSET_MASK or TAO_SUBSET_MATRIXFREE)
16: Output Parameters:
17: . vreduced - the subvector
19: Note:
20: maskvalue should usually be 0.0, unless a pointwise divide will be used.
21: @*/
22: PetscErrorCode VecGetSubVec(Vec vfull, IS is, PetscInt reduced_type, PetscReal maskvalue, Vec *vreduced)
23: {
25: PetscInt nfull,nreduced,nreduced_local,rlow,rhigh,flow,fhigh;
26: PetscInt i,nlocal;
27: PetscReal *fv,*rv;
28: const PetscInt *s;
29: IS ident;
30: VecType vtype;
31: VecScatter scatter;
32: MPI_Comm comm;
38: VecGetSize(vfull, &nfull);
39: ISGetSize(is, &nreduced);
41: if (nreduced == nfull) {
42: VecDestroy(vreduced);
43: VecDuplicate(vfull,vreduced);
44: VecCopy(vfull,*vreduced);
45: } else {
46: switch (reduced_type) {
47: case TAO_SUBSET_SUBVEC:
48: VecGetType(vfull,&vtype);
49: VecGetOwnershipRange(vfull,&flow,&fhigh);
50: ISGetLocalSize(is,&nreduced_local);
51: PetscObjectGetComm((PetscObject)vfull,&comm);
52: if (*vreduced) {
53: VecDestroy(vreduced);
54: }
55: VecCreate(comm,vreduced);
56: VecSetType(*vreduced,vtype);
58: VecSetSizes(*vreduced,nreduced_local,nreduced);
59: VecGetOwnershipRange(*vreduced,&rlow,&rhigh);
60: ISCreateStride(comm,nreduced_local,rlow,1,&ident);
61: VecScatterCreate(vfull,is,*vreduced,ident,&scatter);
62: VecScatterBegin(scatter,vfull,*vreduced,INSERT_VALUES,SCATTER_FORWARD);
63: VecScatterEnd(scatter,vfull,*vreduced,INSERT_VALUES,SCATTER_FORWARD);
64: VecScatterDestroy(&scatter);
65: ISDestroy(&ident);
66: break;
68: case TAO_SUBSET_MASK:
69: case TAO_SUBSET_MATRIXFREE:
70: /* vr[i] = vf[i] if i in is
71: vr[i] = 0 otherwise */
72: if (*vreduced == NULL) {
73: VecDuplicate(vfull,vreduced);
74: }
76: VecSet(*vreduced,maskvalue);
77: ISGetLocalSize(is,&nlocal);
78: VecGetOwnershipRange(vfull,&flow,&fhigh);
79: VecGetArray(vfull,&fv);
80: VecGetArray(*vreduced,&rv);
81: ISGetIndices(is,&s);
82: if (nlocal > (fhigh-flow)) SETERRQ2(PETSC_COMM_WORLD,1,"IS local size %d > Vec local size %d",nlocal,fhigh-flow);
83: for (i=0;i<nlocal;i++) {
84: rv[s[i]-flow] = fv[s[i]-flow];
85: }
86: ISRestoreIndices(is,&s);
87: VecRestoreArray(vfull,&fv);
88: VecRestoreArray(*vreduced,&rv);
89: break;
90: }
91: }
92: return(0);
93: }
97: /*@
98: MatGetSubMat - Gets a submatrix using the IS
100: Input Parameters:
101: + M - the full matrix (n x n)
102: . is - the index set for the submatrix (both row and column index sets need to be the same)
103: . v1 - work vector of dimension n, needed for TAO_SUBSET_MASK option
104: - subset_type - the method TAO is using for subsetting (TAO_SUBSET_SUBVEC, TAO_SUBSET_MASK,
105: TAO_SUBSET_MATRIXFREE)
107: Output Parameters:
108: . Msub - the submatrix
109: @*/
110: PetscErrorCode MatGetSubMat(Mat M, IS is, Vec v1, TaoSubsetType subset_type, Mat *Msub)
111: {
113: IS iscomp;
114: PetscBool flg;
119: MatDestroy(Msub);
120: switch (subset_type) {
121: case TAO_SUBSET_SUBVEC:
122: MatGetSubMatrix(M, is, is, MAT_INITIAL_MATRIX, Msub);
123: break;
125: case TAO_SUBSET_MASK:
126: /* Get Reduced Hessian
127: Msub[i,j] = M[i,j] if i,j in Free_Local or i==j
128: Msub[i,j] = 0 if i!=j and i or j not in Free_Local
129: */
130: PetscOptionsBool("-different_submatrix","use separate hessian matrix when computing submatrices","TaoSubsetType",PETSC_FALSE,&flg,NULL);
131: if (flg == PETSC_TRUE) {
132: MatDuplicate(M, MAT_COPY_VALUES, Msub);
133: } else {
134: /* Act on hessian directly (default) */
135: PetscObjectReference((PetscObject)M);
136: *Msub = M;
137: }
138: /* Save the diagonal to temporary vector */
139: MatGetDiagonal(*Msub,v1);
141: /* Zero out rows and columns */
142: ISComplementVec(is,v1,&iscomp);
144: /* Use v1 instead of 0 here because of PETSc bug */
145: MatZeroRowsColumnsIS(*Msub,iscomp,1.0,v1,v1);
147: ISDestroy(&iscomp);
148: break;
149: case TAO_SUBSET_MATRIXFREE:
150: ISComplementVec(is,v1,&iscomp);
151: MatCreateSubMatrixFree(M,iscomp,iscomp,Msub);
152: ISDestroy(&iscomp);
153: break;
154: }
155: return(0);
156: }