xref: /petsc/src/ksp/pc/impls/deflation/deflationspace.c (revision ccb4e88a40f0b86eaeca07ff64c64e4de2fae686)
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   PetscErrorCode ierr;
47 
48   PetscFunctionBegin;
49   ierr = PetscMalloc1(ncoeffs,&Iidx);CHKERRQ(ierr);
50 
51   ierr = MatCreate(comm,&defl);CHKERRQ(ierr);
52   ierr = MatSetSizes(defl,m,n,M,N);CHKERRQ(ierr);
53   ierr = MatSetUp(defl);CHKERRQ(ierr);
54   ierr = MatSeqAIJSetPreallocation(defl,ncoeffs,NULL);CHKERRQ(ierr);
55   ierr = MatMPIAIJSetPreallocation(defl,ncoeffs,NULL,ncoeffs,NULL);CHKERRQ(ierr);
56   ierr = MatSetOption(defl,MAT_NEW_NONZERO_ALLOCATION_ERR,PETSC_TRUE);CHKERRQ(ierr);
57   ierr = MatSetOption(defl,MAT_NO_OFF_PROC_ENTRIES,PETSC_TRUE);CHKERRQ(ierr);
58 
59   /* Alg 735 Taswell: fvecmat */
60   k = ncoeffs -2;
61   if (trunc) k = k/2;
62 
63   ierr = MatGetOwnershipRange(defl,&ilo,&ihi);CHKERRQ(ierr);
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     ierr = MatSetValues(defl,1,&i,ncoeffs,Iidx,coeffs,INSERT_VALUES);CHKERRQ(ierr);
70     for (j=0; j<ncoeffs; j++) {
71       Iidx[j] += 2;
72       if (Iidx[j] >= N) Iidx[j] = PETSC_MIN_INT;
73     }
74   }
75 
76   ierr = MatAssemblyBegin(defl,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
77   ierr = MatAssemblyEnd(defl,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
78 
79   ierr = PetscFree(Iidx);CHKERRQ(ierr);
80   *H = defl;
81   PetscFunctionReturn(0);
82 }
83 
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;
89   PetscErrorCode ierr;
90 
91   PetscFunctionBegin;
92   /* Haar basis wavelet, level=size */
93   len = pow(2,size);
94   ierr = PetscMalloc2(len,&col,len,&Iidx);CHKERRQ(ierr);
95   val = 1./pow(2,size/2.);
96   for (i=0; i<len; i++) col[i] = val;
97 
98   ierr = PCGetOperators(pc,NULL,&A);CHKERRQ(ierr);
99   ierr = MatGetLocalSize(A,&m,NULL);CHKERRQ(ierr);
100   ierr = MatGetSize(A,&M,NULL);CHKERRQ(ierr);
101   ierr = MatCreate(PetscObjectComm((PetscObject)A),&defl);CHKERRQ(ierr);
102   ierr = MatSetSizes(defl,m,PETSC_DECIDE,M,(PetscInt)ceil(M/(float)len));CHKERRQ(ierr);
103   ierr = MatSetUp(defl);CHKERRQ(ierr);
104   ierr = MatSeqAIJSetPreallocation(defl,size,NULL);CHKERRQ(ierr);
105   ierr = MatMPIAIJSetPreallocation(defl,size,NULL,size,NULL);CHKERRQ(ierr);
106   ierr = MatSetOption(defl,MAT_NEW_NONZERO_ALLOCATION_ERR,PETSC_TRUE);CHKERRQ(ierr);
107 
108   ierr = MatGetOwnershipRangeColumn(defl,&ilo,&ihi);CHKERRQ(ierr);
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     ierr = MatSetValues(defl,len,Iidx,1,&i,col,INSERT_VALUES);CHKERRQ(ierr);
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     ierr = MatSetValues(defl,len,Iidx,1,&ihi,col,INSERT_VALUES);CHKERRQ(ierr);
120   }
121 
122   ierr = MatAssemblyBegin(defl,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
123   ierr = MatAssemblyEnd(defl,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
124 
125   ierr = PetscFree2(col,Iidx);CHKERRQ(ierr);
126   *W = defl;
127   PetscFunctionReturn(0);
128 }
129 
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;
135   PetscErrorCode ierr;
136 
137   PetscFunctionBegin;
138   ierr = PetscObjectGetComm((PetscObject)pc,&comm);CHKERRQ(ierr);
139   ierr = PetscMalloc1(size,&H);CHKERRQ(ierr);
140   ierr = PCGetOperators(pc,&A,NULL);CHKERRQ(ierr);
141   ierr = MatGetLocalSize(A,&m,NULL);CHKERRQ(ierr);
142   ierr = MatGetSize(A,&M,NULL);CHKERRQ(ierr);
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     ierr = PCDeflationCreateSpaceWave(comm,PETSC_DECIDE,m,Mdefl,Ndefl,ncoeffs,coeffs,trunc,&H[i]);CHKERRQ(ierr);
151     ierr = MatGetLocalSize(H[i],&m,NULL);CHKERRQ(ierr);
152     Ndefl = Mdefl;
153   }
154   ierr = MatCreateComposite(comm,size,H,&defl);CHKERRQ(ierr);
155   ierr = MatCompositeSetType(defl,MAT_COMPOSITE_MULTIPLICATIVE);CHKERRQ(ierr);
156   *W = defl;
157 
158   for (i=0; i<size; i++) {
159     ierr = MatDestroy(&H[i]);CHKERRQ(ierr);
160   }
161   ierr = PetscFree(H);CHKERRQ(ierr);
162   PetscFunctionReturn(0);
163 }
164 
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;
172   PetscErrorCode ierr;
173 
174   PetscFunctionBegin;
175   ierr = PCGetOperators(pc,&A,NULL);CHKERRQ(ierr);
176   ierr = MatGetOwnershipRangeColumn(A,&ilo,&ihi);CHKERRQ(ierr);
177   ierr = MatGetSize(A,&M,NULL);CHKERRQ(ierr);
178   ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr);
179   ierr = MPI_Comm_size(comm,&m);CHKERRMPI(ierr);
180   ierr = MatCreate(comm,&defl);CHKERRQ(ierr);
181   ierr = MatSetSizes(defl,ihi-ilo,1,M,m);CHKERRQ(ierr);
182   ierr = MatSetUp(defl);CHKERRQ(ierr);
183   ierr = MatSeqAIJSetPreallocation(defl,1,NULL);CHKERRQ(ierr);
184   ierr = MatMPIAIJSetPreallocation(defl,1,NULL,0,NULL);CHKERRQ(ierr);
185   ierr = MatSetOption(defl,MAT_NEW_NONZERO_ALLOCATION_ERR,PETSC_TRUE);CHKERRQ(ierr);
186   ierr = MatSetOption(defl,MAT_NO_OFF_PROC_ENTRIES,PETSC_TRUE);CHKERRQ(ierr);
187 
188   ierr = PetscMalloc2(ihi-ilo,&col,ihi-ilo,&Iidx);CHKERRQ(ierr);
189   for (i=ilo; i<ihi; i++) {
190     Iidx[i-ilo] = i;
191     col[i-ilo] = 1;
192   }
193   ierr = MPI_Comm_rank(comm,&m);CHKERRMPI(ierr);
194   i = m;
195   ierr = MatSetValues(defl,ihi-ilo,Iidx,1,&i,col,INSERT_VALUES);CHKERRQ(ierr);
196 
197   ierr = MatAssemblyBegin(defl,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
198   ierr = MatAssemblyEnd(defl,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
199 
200   ierr = PetscFree2(col,Iidx);CHKERRQ(ierr);
201   *W = defl;
202   PetscFunctionReturn(0);
203 }
204 
205 PetscErrorCode PCDeflationComputeSpace(PC pc)
206 {
207   Mat            defl;
208   PetscBool      transp=PETSC_TRUE;
209   PC_Deflation   *def = (PC_Deflation*)pc->data;
210   PetscErrorCode ierr;
211 
212   PetscFunctionBegin;
213   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
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       ierr = PCDeflationGetSpaceHaar(pc,&defl,def->spacesize);CHKERRQ(ierr);break;
219     case PC_DEFLATION_SPACE_DB2:
220       ierr = PCDeflationGetSpaceWave(pc,&defl,def->spacesize,2,db2,PetscNot(def->extendsp));CHKERRQ(ierr);break;
221     case PC_DEFLATION_SPACE_DB4:
222       ierr = PCDeflationGetSpaceWave(pc,&defl,def->spacesize,4,db4,PetscNot(def->extendsp));CHKERRQ(ierr);break;
223     case PC_DEFLATION_SPACE_DB8:
224       ierr = PCDeflationGetSpaceWave(pc,&defl,def->spacesize,8,db8,PetscNot(def->extendsp));CHKERRQ(ierr);break;
225     case PC_DEFLATION_SPACE_DB16:
226       ierr = PCDeflationGetSpaceWave(pc,&defl,def->spacesize,16,db16,PetscNot(def->extendsp));CHKERRQ(ierr);break;
227     case PC_DEFLATION_SPACE_BIORTH22:
228       ierr = PCDeflationGetSpaceWave(pc,&defl,def->spacesize,6,biorth22,PetscNot(def->extendsp));CHKERRQ(ierr);break;
229     case PC_DEFLATION_SPACE_MEYER:
230       ierr = PCDeflationGetSpaceWave(pc,&defl,def->spacesize,62,meyer,PetscNot(def->extendsp));CHKERRQ(ierr);break;
231     case PC_DEFLATION_SPACE_AGGREGATION:
232       transp = PETSC_FALSE;
233       ierr = PCDeflationGetSpaceAggregation(pc,&defl);CHKERRQ(ierr);break;
234     default: SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONG,"Wrong PCDeflationSpaceType specified");
235   }
236 
237   ierr = PCDeflationSetSpace(pc,defl,transp);CHKERRQ(ierr);
238   ierr = MatDestroy(&defl);CHKERRQ(ierr);
239   PetscFunctionReturn(0);
240 }
241 
242