xref: /petsc/src/ksp/pc/impls/deflation/deflationspace.c (revision dfd676b1a855b7f967ece75a22ee7f6626d10f89)
1 #include <../src/ksp/pc/impls/deflation/deflation.h> /*I "petscksp.h" I*/
2 
3 PetscScalar db2[] = {0.7071067811865476,0.7071067811865476};
4 
5 PetscScalar db4[] = {-0.12940952255092145,0.22414386804185735,0.836516303737469,0.48296291314469025};
6 
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};
15 
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};
32 
33 PetscScalar biorth22[] = {0.0,
34 -0.1767766952966369,
35 0.3535533905932738,
36 1.0606601717798214,
37 0.3535533905932738,
38 -0.1767766952966369};
39 
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};
41 
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;
46 
47   PetscFunctionBegin;
48   PetscCall(PetscMalloc1(ncoeffs,&Iidx));
49 
50   PetscCall(MatCreate(comm,&defl));
51   PetscCall(MatSetSizes(defl,m,n,M,N));
52   PetscCall(MatSetUp(defl));
53   PetscCall(MatSeqAIJSetPreallocation(defl,ncoeffs,NULL));
54   PetscCall(MatMPIAIJSetPreallocation(defl,ncoeffs,NULL,ncoeffs,NULL));
55   PetscCall(MatSetOption(defl,MAT_NEW_NONZERO_ALLOCATION_ERR,PETSC_TRUE));
56   PetscCall(MatSetOption(defl,MAT_NO_OFF_PROC_ENTRIES,PETSC_TRUE));
57 
58   /* Alg 735 Taswell: fvecmat */
59   k = ncoeffs -2;
60   if (trunc) k = k/2;
61 
62   PetscCall(MatGetOwnershipRange(defl,&ilo,&ihi));
63   for (i=0; i<ncoeffs; i++) {
64     Iidx[i] = i+ilo*2 -k;
65     if (Iidx[i] >= N) Iidx[i] = PETSC_MIN_INT;
66   }
67   for (i=ilo; i<ihi; i++) {
68     PetscCall(MatSetValues(defl,1,&i,ncoeffs,Iidx,coeffs,INSERT_VALUES));
69     for (j=0; j<ncoeffs; j++) {
70       Iidx[j] += 2;
71       if (Iidx[j] >= N) Iidx[j] = PETSC_MIN_INT;
72     }
73   }
74 
75   PetscCall(MatAssemblyBegin(defl,MAT_FINAL_ASSEMBLY));
76   PetscCall(MatAssemblyEnd(defl,MAT_FINAL_ASSEMBLY));
77 
78   PetscCall(PetscFree(Iidx));
79   *H = defl;
80   PetscFunctionReturn(0);
81 }
82 
83 PetscErrorCode PCDeflationGetSpaceHaar(PC pc,Mat *W,PetscInt size)
84 {
85   Mat            A,defl;
86   PetscInt       i,j,len,ilo,ihi,*Iidx,m,M;
87   PetscScalar    *col,val;
88 
89   PetscFunctionBegin;
90   /* Haar basis wavelet, level=size */
91   len = pow(2,size);
92   PetscCall(PetscMalloc2(len,&col,len,&Iidx));
93   val = 1./pow(2,size/2.);
94   for (i=0; i<len; i++) col[i] = val;
95 
96   PetscCall(PCGetOperators(pc,NULL,&A));
97   PetscCall(MatGetLocalSize(A,&m,NULL));
98   PetscCall(MatGetSize(A,&M,NULL));
99   PetscCall(MatCreate(PetscObjectComm((PetscObject)A),&defl));
100   PetscCall(MatSetSizes(defl,m,PETSC_DECIDE,M,PetscCeilInt(M,len)));
101   PetscCall(MatSetUp(defl));
102   PetscCall(MatSeqAIJSetPreallocation(defl,size,NULL));
103   PetscCall(MatMPIAIJSetPreallocation(defl,size,NULL,size,NULL));
104   PetscCall(MatSetOption(defl,MAT_NEW_NONZERO_ALLOCATION_ERR,PETSC_TRUE));
105 
106   PetscCall(MatGetOwnershipRangeColumn(defl,&ilo,&ihi));
107   for (i=0; i<len; i++) Iidx[i] = i+ilo*len;
108   if (M%len && ihi == PetscCeilInt(M,len)) ihi -= 1;
109   for (i=ilo; i<ihi; i++) {
110     PetscCall(MatSetValues(defl,len,Iidx,1,&i,col,INSERT_VALUES));
111     for (j=0; j<len; j++) Iidx[j] += len;
112   }
113   if (M%len && ihi+1 == PetscCeilInt(M,len)) {
114     len = M%len;
115     val = 1./pow(pow(2,len),0.5);
116     for (i=0; i<len; i++) col[i] = val;
117     PetscCall(MatSetValues(defl,len,Iidx,1,&ihi,col,INSERT_VALUES));
118   }
119 
120   PetscCall(MatAssemblyBegin(defl,MAT_FINAL_ASSEMBLY));
121   PetscCall(MatAssemblyEnd(defl,MAT_FINAL_ASSEMBLY));
122 
123   PetscCall(PetscFree2(col,Iidx));
124   *W = defl;
125   PetscFunctionReturn(0);
126 }
127 
128 PetscErrorCode PCDeflationGetSpaceWave(PC pc,Mat *W,PetscInt size,PetscInt ncoeffs,PetscScalar *coeffs,PetscBool trunc)
129 {
130   Mat            A,*H,defl;
131   PetscInt       i,m,M,Mdefl,Ndefl;
132   MPI_Comm       comm;
133 
134   PetscFunctionBegin;
135   PetscCall(PetscObjectGetComm((PetscObject)pc,&comm));
136   PetscCall(PetscMalloc1(size,&H));
137   PetscCall(PCGetOperators(pc,&A,NULL));
138   PetscCall(MatGetLocalSize(A,&m,NULL));
139   PetscCall(MatGetSize(A,&M,NULL));
140   Mdefl = M;
141   Ndefl = M;
142   for (i=0; i<size; i++) {
143     if (Mdefl%2)  {
144       if (trunc) Mdefl = (PetscInt)PetscCeilReal(Mdefl/2.);
145       else       Mdefl = (PetscInt)PetscFloorReal((ncoeffs+Mdefl-1)/2.);
146     } else       Mdefl = Mdefl/2;
147     PetscCall(PCDeflationCreateSpaceWave(comm,PETSC_DECIDE,m,Mdefl,Ndefl,ncoeffs,coeffs,trunc,&H[i]));
148     PetscCall(MatGetLocalSize(H[i],&m,NULL));
149     Ndefl = Mdefl;
150   }
151   PetscCall(MatCreateComposite(comm,size,H,&defl));
152   PetscCall(MatCompositeSetType(defl,MAT_COMPOSITE_MULTIPLICATIVE));
153   *W = defl;
154 
155   for (i=0; i<size; i++) {
156     PetscCall(MatDestroy(&H[i]));
157   }
158   PetscCall(PetscFree(H));
159   PetscFunctionReturn(0);
160 }
161 
162 PetscErrorCode PCDeflationGetSpaceAggregation(PC pc,Mat *W)
163 {
164   Mat            A,defl;
165   PetscInt       i,ilo,ihi,*Iidx,M;
166   PetscMPIInt    m;
167   PetscScalar    *col;
168   MPI_Comm       comm;
169 
170   PetscFunctionBegin;
171   PetscCall(PCGetOperators(pc,&A,NULL));
172   PetscCall(MatGetOwnershipRangeColumn(A,&ilo,&ihi));
173   PetscCall(MatGetSize(A,&M,NULL));
174   PetscCall(PetscObjectGetComm((PetscObject)A,&comm));
175   PetscCallMPI(MPI_Comm_size(comm,&m));
176   PetscCall(MatCreate(comm,&defl));
177   PetscCall(MatSetSizes(defl,ihi-ilo,1,M,m));
178   PetscCall(MatSetUp(defl));
179   PetscCall(MatSeqAIJSetPreallocation(defl,1,NULL));
180   PetscCall(MatMPIAIJSetPreallocation(defl,1,NULL,0,NULL));
181   PetscCall(MatSetOption(defl,MAT_NEW_NONZERO_ALLOCATION_ERR,PETSC_TRUE));
182   PetscCall(MatSetOption(defl,MAT_NO_OFF_PROC_ENTRIES,PETSC_TRUE));
183 
184   PetscCall(PetscMalloc2(ihi-ilo,&col,ihi-ilo,&Iidx));
185   for (i=ilo; i<ihi; i++) {
186     Iidx[i-ilo] = i;
187     col[i-ilo] = 1;
188   }
189   PetscCallMPI(MPI_Comm_rank(comm,&m));
190   i = m;
191   PetscCall(MatSetValues(defl,ihi-ilo,Iidx,1,&i,col,INSERT_VALUES));
192 
193   PetscCall(MatAssemblyBegin(defl,MAT_FINAL_ASSEMBLY));
194   PetscCall(MatAssemblyEnd(defl,MAT_FINAL_ASSEMBLY));
195 
196   PetscCall(PetscFree2(col,Iidx));
197   *W = defl;
198   PetscFunctionReturn(0);
199 }
200 
201 PetscErrorCode PCDeflationComputeSpace(PC pc)
202 {
203   Mat            defl;
204   PetscBool      transp=PETSC_TRUE;
205   PC_Deflation   *def = (PC_Deflation*)pc->data;
206 
207   PetscFunctionBegin;
208   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
209   PetscCheck(def->spacesize >= 1,PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONG,"Wrong PCDeflation space size specified: %" PetscInt_FMT,def->spacesize);
210   switch (def->spacetype) {
211     case PC_DEFLATION_SPACE_HAAR:
212       transp = PETSC_FALSE;
213       PetscCall(PCDeflationGetSpaceHaar(pc,&defl,def->spacesize));break;
214     case PC_DEFLATION_SPACE_DB2:
215       PetscCall(PCDeflationGetSpaceWave(pc,&defl,def->spacesize,2,db2,PetscNot(def->extendsp)));break;
216     case PC_DEFLATION_SPACE_DB4:
217       PetscCall(PCDeflationGetSpaceWave(pc,&defl,def->spacesize,4,db4,PetscNot(def->extendsp)));break;
218     case PC_DEFLATION_SPACE_DB8:
219       PetscCall(PCDeflationGetSpaceWave(pc,&defl,def->spacesize,8,db8,PetscNot(def->extendsp)));break;
220     case PC_DEFLATION_SPACE_DB16:
221       PetscCall(PCDeflationGetSpaceWave(pc,&defl,def->spacesize,16,db16,PetscNot(def->extendsp)));break;
222     case PC_DEFLATION_SPACE_BIORTH22:
223       PetscCall(PCDeflationGetSpaceWave(pc,&defl,def->spacesize,6,biorth22,PetscNot(def->extendsp)));break;
224     case PC_DEFLATION_SPACE_MEYER:
225       PetscCall(PCDeflationGetSpaceWave(pc,&defl,def->spacesize,62,meyer,PetscNot(def->extendsp)));break;
226     case PC_DEFLATION_SPACE_AGGREGATION:
227       transp = PETSC_FALSE;
228       PetscCall(PCDeflationGetSpaceAggregation(pc,&defl));break;
229     default: SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONG,"Wrong PCDeflationSpaceType specified");
230   }
231 
232   PetscCall(PCDeflationSetSpace(pc,defl,transp));
233   PetscCall(MatDestroy(&defl));
234   PetscFunctionReturn(0);
235 }
236