Actual source code: deflationspace.c

petsc-3.14.3 2021-01-09
Report Typos and Errors
  1: #include <../src/ksp/pc/impls/deflation/deflation.h>

  3: PetscScalar db2[] = {0.7071067811865476,0.7071067811865476};

  5: PetscScalar db4[] = {-0.12940952255092145,0.22414386804185735,0.836516303737469,0.48296291314469025};

  7: PetscScalar db8[] = {-0.010597401784997278,
  8: 0.032883011666982945,
  9: 0.030841381835986965,
 10: -0.18703481171888114,
 11: -0.02798376941698385,
 12: 0.6308807679295904,
 13: 0.7148465705525415,
 14: 0.23037781330885523};

 16: PetscScalar db16[] = {-0.00011747678400228192,
 17: 0.0006754494059985568,
 18: -0.0003917403729959771,
 19: -0.00487035299301066,
 20: 0.008746094047015655,
 21: 0.013981027917015516,
 22: -0.04408825393106472,
 23: -0.01736930100202211,
 24: 0.128747426620186,
 25: 0.00047248457399797254,
 26: -0.2840155429624281,
 27: -0.015829105256023893,
 28: 0.5853546836548691,
 29: 0.6756307362980128,
 30: 0.3128715909144659,
 31: 0.05441584224308161};

 33: PetscScalar biorth22[] = {0.0,
 34: -0.1767766952966369,
 35: 0.3535533905932738,
 36: 1.0606601717798214,
 37: 0.3535533905932738,
 38: -0.1767766952966369};

 40: PetscScalar meyer[] = {0.0,-1.009999956941423e-12,8.519459636796214e-09,-1.111944952595278e-08,-1.0798819539621958e-08,6.066975741351135e-08,-1.0866516536735883e-07,8.200680650386481e-08,1.1783004497663934e-07,-5.506340565252278e-07,1.1307947017916706e-06,-1.489549216497156e-06,7.367572885903746e-07,3.20544191334478e-06,-1.6312699734552807e-05,6.554305930575149e-05,-0.0006011502343516092,-0.002704672124643725,0.002202534100911002,0.006045814097323304,-0.006387718318497156,-0.011061496392513451,0.015270015130934803,0.017423434103729693,-0.03213079399021176,-0.024348745906078023,0.0637390243228016,0.030655091960824263,-0.13284520043622938,-0.035087555656258346,0.44459300275757724,0.7445855923188063,0.44459300275757724,-0.035087555656258346,-0.13284520043622938,0.030655091960824263,0.0637390243228016,-0.024348745906078023,-0.03213079399021176,0.017423434103729693,0.015270015130934803,-0.011061496392513451,-0.006387718318497156,0.006045814097323304,0.002202534100911002,-0.002704672124643725,-0.0006011502343516092,6.554305930575149e-05,-1.6312699734552807e-05,3.20544191334478e-06,7.367572885903746e-07,-1.489549216497156e-06,1.1307947017916706e-06,-5.506340565252278e-07,1.1783004497663934e-07,8.200680650386481e-08,-1.0866516536735883e-07,6.066975741351135e-08,-1.0798819539621958e-08,-1.111944952595278e-08,8.519459636796214e-09,-1.009999956941423e-12};

 42: static PetscErrorCode PCDeflationCreateSpaceWave(MPI_Comm comm,PetscInt m,PetscInt n,PetscInt M,PetscInt N,PetscInt ncoeffs,PetscScalar *coeffs,PetscBool trunc,Mat *H)
 43: {
 44:   Mat            defl;
 45:   PetscInt       i,j,k,ilo,ihi,*Iidx;

 49:   PetscMalloc1(ncoeffs,&Iidx);

 51:   MatCreate(comm,&defl);
 52:   MatSetSizes(defl,m,n,M,N);
 53:   MatSetUp(defl);
 54:   MatSeqAIJSetPreallocation(defl,ncoeffs,NULL);
 55:   MatMPIAIJSetPreallocation(defl,ncoeffs,NULL,ncoeffs,NULL);
 56:   MatSetOption(defl,MAT_NEW_NONZERO_ALLOCATION_ERR,PETSC_TRUE);
 57:   MatSetOption(defl,MAT_NO_OFF_PROC_ENTRIES,PETSC_TRUE);

 59:   /* Alg 735 Taswell: fvecmat */
 60:   k = ncoeffs -2;
 61:   if (trunc) k = k/2;

 63:   MatGetOwnershipRange(defl,&ilo,&ihi);
 64:   for (i=0; i<ncoeffs; i++) {
 65:     Iidx[i] = i+ilo*2 -k;
 66:     if (Iidx[i] >= N) Iidx[i] = PETSC_MIN_INT;
 67:   }
 68:   for (i=ilo; i<ihi; i++) {
 69:     MatSetValues(defl,1,&i,ncoeffs,Iidx,coeffs,INSERT_VALUES);
 70:     for (j=0; j<ncoeffs; j++) {
 71:       Iidx[j] += 2;
 72:       if (Iidx[j] >= N) Iidx[j] = PETSC_MIN_INT;
 73:     }
 74:   }

 76:   MatAssemblyBegin(defl,MAT_FINAL_ASSEMBLY);
 77:   MatAssemblyEnd(defl,MAT_FINAL_ASSEMBLY);

 79:   PetscFree(Iidx);
 80:   *H = defl;
 81:   return(0);
 82: }

 84: PetscErrorCode PCDeflationGetSpaceHaar(PC pc,Mat *W,PetscInt size)
 85: {
 86:   Mat            A,defl;
 87:   PetscInt       i,j,len,ilo,ihi,*Iidx,m,M;
 88:   PetscScalar    *col,val;

 92:   /* Haar basis wavelet, level=size */
 93:   len = pow(2,size);
 94:   PetscMalloc2(len,&col,len,&Iidx);
 95:   val = 1./pow(2,size/2.);
 96:   for (i=0; i<len; i++) col[i] = val;

 98:   PCGetOperators(pc,NULL,&A);
 99:   MatGetLocalSize(A,&m,NULL);
100:   MatGetSize(A,&M,NULL);
101:   MatCreate(PetscObjectComm((PetscObject)A),&defl);
102:   MatSetSizes(defl,m,PETSC_DECIDE,M,(PetscInt)ceil(M/(float)len));
103:   MatSetUp(defl);
104:   MatSeqAIJSetPreallocation(defl,size,NULL);
105:   MatMPIAIJSetPreallocation(defl,size,NULL,size,NULL);
106:   MatSetOption(defl,MAT_NEW_NONZERO_ALLOCATION_ERR,PETSC_TRUE);

108:   MatGetOwnershipRangeColumn(defl,&ilo,&ihi);
109:   for (i=0; i<len; i++) Iidx[i] = i+ilo*len;
110:   if (M%len && ihi == (int)ceil(M/(float)len)) ihi -= 1;
111:   for (i=ilo; i<ihi; i++) {
112:     MatSetValues(defl,len,Iidx,1,&i,col,INSERT_VALUES);
113:     for (j=0; j<len; j++) Iidx[j] += len;
114:   }
115:   if (M%len && ihi+1 == ceil(M/(float)len)) {
116:     len = M%len;
117:     val = 1./pow(pow(2,len),0.5);
118:     for (i=0; i<len; i++) col[i] = val;
119:     MatSetValues(defl,len,Iidx,1,&ihi,col,INSERT_VALUES);
120:   }

122:   MatAssemblyBegin(defl,MAT_FINAL_ASSEMBLY);
123:   MatAssemblyEnd(defl,MAT_FINAL_ASSEMBLY);

125:   PetscFree2(col,Iidx);
126:   *W = defl;
127:   return(0);
128: }

130: PetscErrorCode PCDeflationGetSpaceWave(PC pc,Mat *W,PetscInt size,PetscInt ncoeffs,PetscScalar *coeffs,PetscBool trunc)
131: {
132:   Mat            A,*H,defl;
133:   PetscInt       i,m,M,Mdefl,Ndefl;
134:   MPI_Comm       comm;

138:   PetscObjectGetComm((PetscObject)pc,&comm);
139:   PetscMalloc1(size,&H);
140:   PCGetOperators(pc,&A,NULL);
141:   MatGetLocalSize(A,&m,NULL);
142:   MatGetSize(A,&M,NULL);
143:   Mdefl = M;
144:   Ndefl = M;
145:   for (i=0; i<size; i++) {
146:     if (Mdefl%2)  {
147:       if (trunc) Mdefl = (PetscInt)PetscCeilReal(Mdefl/2.);
148:       else       Mdefl = (PetscInt)PetscFloorReal((ncoeffs+Mdefl-1)/2.);
149:     } else       Mdefl = Mdefl/2;
150:     PCDeflationCreateSpaceWave(comm,PETSC_DECIDE,m,Mdefl,Ndefl,ncoeffs,coeffs,trunc,&H[i]);
151:     MatGetLocalSize(H[i],&m,NULL);
152:     Ndefl = Mdefl;
153:   }
154:   MatCreateComposite(comm,size,H,&defl);
155:   MatCompositeSetType(defl,MAT_COMPOSITE_MULTIPLICATIVE);
156:   *W = defl;

158:   for (i=0; i<size; i++) {
159:     MatDestroy(&H[i]);
160:   }
161:   PetscFree(H);
162:   return(0);
163: }

165: PetscErrorCode PCDeflationGetSpaceAggregation(PC pc,Mat *W)
166: {
167:   Mat            A,defl;
168:   PetscInt       i,ilo,ihi,*Iidx,M;
169:   PetscMPIInt    m;
170:   PetscScalar    *col;
171:   MPI_Comm       comm;

175:   PCGetOperators(pc,&A,NULL);
176:   MatGetOwnershipRangeColumn(A,&ilo,&ihi);
177:   MatGetSize(A,&M,NULL);
178:   PetscObjectGetComm((PetscObject)A,&comm);
179:   MPI_Comm_size(comm,&m);
180:   MatCreate(comm,&defl);
181:   MatSetSizes(defl,ihi-ilo,1,M,m);
182:   MatSetUp(defl);
183:   MatSeqAIJSetPreallocation(defl,1,NULL);
184:   MatMPIAIJSetPreallocation(defl,1,NULL,0,NULL);
185:   MatSetOption(defl,MAT_NEW_NONZERO_ALLOCATION_ERR,PETSC_TRUE);
186:   MatSetOption(defl,MAT_NO_OFF_PROC_ENTRIES,PETSC_TRUE);

188:   PetscMalloc2(ihi-ilo,&col,ihi-ilo,&Iidx);
189:   for (i=ilo; i<ihi; i++) {
190:     Iidx[i-ilo] = i;
191:     col[i-ilo] = 1;
192:   }
193:   MPI_Comm_rank(comm,&m);
194:   i = m;
195:   MatSetValues(defl,ihi-ilo,Iidx,1,&i,col,INSERT_VALUES);

197:   MatAssemblyBegin(defl,MAT_FINAL_ASSEMBLY);
198:   MatAssemblyEnd(defl,MAT_FINAL_ASSEMBLY);

200:   PetscFree2(col,Iidx);
201:   *W = defl;
202:   return(0);
203: }

205: PetscErrorCode PCDeflationComputeSpace(PC pc)
206: {
207:   Mat            defl;
208:   PetscBool      transp=PETSC_TRUE;
209:   PC_Deflation   *def = (PC_Deflation*)pc->data;

214:   if (def->spacesize < 1) SETERRQ1(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONG,"Wrong PCDeflation space size specified: %D",def->spacesize);
215:   switch (def->spacetype) {
216:     case PC_DEFLATION_SPACE_HAAR:
217:       transp = PETSC_FALSE;
218:       PCDeflationGetSpaceHaar(pc,&defl,def->spacesize);break;
219:     case PC_DEFLATION_SPACE_DB2:
220:       PCDeflationGetSpaceWave(pc,&defl,def->spacesize,2,db2,PetscNot(def->extendsp));break;
221:     case PC_DEFLATION_SPACE_DB4:
222:       PCDeflationGetSpaceWave(pc,&defl,def->spacesize,4,db4,PetscNot(def->extendsp));break;
223:     case PC_DEFLATION_SPACE_DB8:
224:       PCDeflationGetSpaceWave(pc,&defl,def->spacesize,8,db8,PetscNot(def->extendsp));break;
225:     case PC_DEFLATION_SPACE_DB16:
226:       PCDeflationGetSpaceWave(pc,&defl,def->spacesize,16,db16,PetscNot(def->extendsp));break;
227:     case PC_DEFLATION_SPACE_BIORTH22:
228:       PCDeflationGetSpaceWave(pc,&defl,def->spacesize,6,biorth22,PetscNot(def->extendsp));break;
229:     case PC_DEFLATION_SPACE_MEYER:
230:       PCDeflationGetSpaceWave(pc,&defl,def->spacesize,62,meyer,PetscNot(def->extendsp));break;
231:     case PC_DEFLATION_SPACE_AGGREGATION:
232:       transp = PETSC_FALSE;
233:       PCDeflationGetSpaceAggregation(pc,&defl);break;
234:     default: SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONG,"Wrong PCDeflationSpaceType specified");
235:   }

237:   PCDeflationSetSpace(pc,defl,transp);
238:   MatDestroy(&defl);
239:   return(0);
240: }