15170378fSJakub Kruzik #include <../src/ksp/pc/impls/deflation/deflation.h> /*I "petscksp.h" I*/
2e53e0a0dSJakub Kruzik
3e53e0a0dSJakub Kruzik PetscScalar db2[] = {0.7071067811865476, 0.7071067811865476};
4e53e0a0dSJakub Kruzik
5e53e0a0dSJakub Kruzik PetscScalar db4[] = {-0.12940952255092145, 0.22414386804185735, 0.836516303737469, 0.48296291314469025};
6e53e0a0dSJakub Kruzik
79371c9d4SSatish Balay PetscScalar db8[] = {-0.010597401784997278, 0.032883011666982945, 0.030841381835986965, -0.18703481171888114, -0.02798376941698385, 0.6308807679295904, 0.7148465705525415, 0.23037781330885523};
8e53e0a0dSJakub Kruzik
99371c9d4SSatish Balay PetscScalar db16[] = {-0.00011747678400228192, 0.0006754494059985568, -0.0003917403729959771, -0.00487035299301066, 0.008746094047015655, 0.013981027917015516, -0.04408825393106472, -0.01736930100202211,
109371c9d4SSatish Balay 0.128747426620186, 0.00047248457399797254, -0.2840155429624281, -0.015829105256023893, 0.5853546836548691, 0.6756307362980128, 0.3128715909144659, 0.05441584224308161};
11e53e0a0dSJakub Kruzik
129371c9d4SSatish Balay PetscScalar biorth22[] = {0.0, -0.1767766952966369, 0.3535533905932738, 1.0606601717798214, 0.3535533905932738, -0.1767766952966369};
13e53e0a0dSJakub Kruzik
1465149469SJakub Kruzik 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};
15e53e0a0dSJakub Kruzik
PCDeflationCreateSpaceWave(MPI_Comm comm,PetscInt m,PetscInt n,PetscInt M,PetscInt N,PetscInt ncoeffs,PetscScalar * coeffs,PetscBool trunc,Mat * H)16d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCDeflationCreateSpaceWave(MPI_Comm comm, PetscInt m, PetscInt n, PetscInt M, PetscInt N, PetscInt ncoeffs, PetscScalar *coeffs, PetscBool trunc, Mat *H)
17d71ae5a4SJacob Faibussowitsch {
18e53e0a0dSJakub Kruzik Mat defl;
19e53e0a0dSJakub Kruzik PetscInt i, j, k, ilo, ihi, *Iidx;
20e53e0a0dSJakub Kruzik
21e53e0a0dSJakub Kruzik PetscFunctionBegin;
229566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(ncoeffs, &Iidx));
23e53e0a0dSJakub Kruzik
249566063dSJacob Faibussowitsch PetscCall(MatCreate(comm, &defl));
259566063dSJacob Faibussowitsch PetscCall(MatSetSizes(defl, m, n, M, N));
269566063dSJacob Faibussowitsch PetscCall(MatSetUp(defl));
279566063dSJacob Faibussowitsch PetscCall(MatSeqAIJSetPreallocation(defl, ncoeffs, NULL));
289566063dSJacob Faibussowitsch PetscCall(MatMPIAIJSetPreallocation(defl, ncoeffs, NULL, ncoeffs, NULL));
299566063dSJacob Faibussowitsch PetscCall(MatSetOption(defl, MAT_NEW_NONZERO_ALLOCATION_ERR, PETSC_TRUE));
309566063dSJacob Faibussowitsch PetscCall(MatSetOption(defl, MAT_NO_OFF_PROC_ENTRIES, PETSC_TRUE));
31e53e0a0dSJakub Kruzik
32e53e0a0dSJakub Kruzik /* Alg 735 Taswell: fvecmat */
33e53e0a0dSJakub Kruzik k = ncoeffs - 2;
34e53e0a0dSJakub Kruzik if (trunc) k = k / 2;
35e53e0a0dSJakub Kruzik
369566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(defl, &ilo, &ihi));
37e53e0a0dSJakub Kruzik for (i = 0; i < ncoeffs; i++) {
38e53e0a0dSJakub Kruzik Iidx[i] = i + ilo * 2 - k;
39*1690c2aeSBarry Smith if (Iidx[i] >= N) Iidx[i] = PETSC_INT_MIN;
40e53e0a0dSJakub Kruzik }
41e53e0a0dSJakub Kruzik for (i = ilo; i < ihi; i++) {
429566063dSJacob Faibussowitsch PetscCall(MatSetValues(defl, 1, &i, ncoeffs, Iidx, coeffs, INSERT_VALUES));
43e53e0a0dSJakub Kruzik for (j = 0; j < ncoeffs; j++) {
44e53e0a0dSJakub Kruzik Iidx[j] += 2;
45*1690c2aeSBarry Smith if (Iidx[j] >= N) Iidx[j] = PETSC_INT_MIN;
46e53e0a0dSJakub Kruzik }
47e53e0a0dSJakub Kruzik }
48e53e0a0dSJakub Kruzik
499566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(defl, MAT_FINAL_ASSEMBLY));
509566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(defl, MAT_FINAL_ASSEMBLY));
51e53e0a0dSJakub Kruzik
529566063dSJacob Faibussowitsch PetscCall(PetscFree(Iidx));
53e53e0a0dSJakub Kruzik *H = defl;
543ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
55e53e0a0dSJakub Kruzik }
56e53e0a0dSJakub Kruzik
PCDeflationGetSpaceHaar(PC pc,Mat * W,PetscInt size)5766976f2fSJacob Faibussowitsch static PetscErrorCode PCDeflationGetSpaceHaar(PC pc, Mat *W, PetscInt size)
58d71ae5a4SJacob Faibussowitsch {
59e53e0a0dSJakub Kruzik Mat A, defl;
60e53e0a0dSJakub Kruzik PetscInt i, j, len, ilo, ihi, *Iidx, m, M;
619e56ec8aSJakub Kruzik PetscScalar *col, val;
62e53e0a0dSJakub Kruzik
63e53e0a0dSJakub Kruzik PetscFunctionBegin;
64e53e0a0dSJakub Kruzik /* Haar basis wavelet, level=size */
65e53e0a0dSJakub Kruzik len = pow(2, size);
669566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(len, &col, len, &Iidx));
67e53e0a0dSJakub Kruzik val = 1. / pow(2, size / 2.);
68e53e0a0dSJakub Kruzik for (i = 0; i < len; i++) col[i] = val;
69e53e0a0dSJakub Kruzik
709566063dSJacob Faibussowitsch PetscCall(PCGetOperators(pc, NULL, &A));
719566063dSJacob Faibussowitsch PetscCall(MatGetLocalSize(A, &m, NULL));
729566063dSJacob Faibussowitsch PetscCall(MatGetSize(A, &M, NULL));
739566063dSJacob Faibussowitsch PetscCall(MatCreate(PetscObjectComm((PetscObject)A), &defl));
749566063dSJacob Faibussowitsch PetscCall(MatSetSizes(defl, m, PETSC_DECIDE, M, PetscCeilInt(M, len)));
759566063dSJacob Faibussowitsch PetscCall(MatSetUp(defl));
769566063dSJacob Faibussowitsch PetscCall(MatSeqAIJSetPreallocation(defl, size, NULL));
779566063dSJacob Faibussowitsch PetscCall(MatMPIAIJSetPreallocation(defl, size, NULL, size, NULL));
789566063dSJacob Faibussowitsch PetscCall(MatSetOption(defl, MAT_NEW_NONZERO_ALLOCATION_ERR, PETSC_TRUE));
79e53e0a0dSJakub Kruzik
809566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRangeColumn(defl, &ilo, &ihi));
81e53e0a0dSJakub Kruzik for (i = 0; i < len; i++) Iidx[i] = i + ilo * len;
82faa75363SBarry Smith if (M % len && ihi == PetscCeilInt(M, len)) ihi -= 1;
83e53e0a0dSJakub Kruzik for (i = ilo; i < ihi; i++) {
849566063dSJacob Faibussowitsch PetscCall(MatSetValues(defl, len, Iidx, 1, &i, col, INSERT_VALUES));
85e53e0a0dSJakub Kruzik for (j = 0; j < len; j++) Iidx[j] += len;
86e53e0a0dSJakub Kruzik }
87faa75363SBarry Smith if (M % len && ihi + 1 == PetscCeilInt(M, len)) {
88e53e0a0dSJakub Kruzik len = M % len;
89e53e0a0dSJakub Kruzik val = 1. / pow(pow(2, len), 0.5);
90e53e0a0dSJakub Kruzik for (i = 0; i < len; i++) col[i] = val;
919566063dSJacob Faibussowitsch PetscCall(MatSetValues(defl, len, Iidx, 1, &ihi, col, INSERT_VALUES));
92e53e0a0dSJakub Kruzik }
93e53e0a0dSJakub Kruzik
949566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(defl, MAT_FINAL_ASSEMBLY));
959566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(defl, MAT_FINAL_ASSEMBLY));
96e53e0a0dSJakub Kruzik
979566063dSJacob Faibussowitsch PetscCall(PetscFree2(col, Iidx));
98e53e0a0dSJakub Kruzik *W = defl;
993ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
100e53e0a0dSJakub Kruzik }
101e53e0a0dSJakub Kruzik
PCDeflationGetSpaceWave(PC pc,Mat * W,PetscInt size,PetscInt ncoeffs,PetscScalar * coeffs,PetscBool trunc)10266976f2fSJacob Faibussowitsch static PetscErrorCode PCDeflationGetSpaceWave(PC pc, Mat *W, PetscInt size, PetscInt ncoeffs, PetscScalar *coeffs, PetscBool trunc)
103d71ae5a4SJacob Faibussowitsch {
104e53e0a0dSJakub Kruzik Mat A, *H, defl;
105e53e0a0dSJakub Kruzik PetscInt i, m, M, Mdefl, Ndefl;
106e53e0a0dSJakub Kruzik MPI_Comm comm;
107e53e0a0dSJakub Kruzik
108e53e0a0dSJakub Kruzik PetscFunctionBegin;
1099566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)pc, &comm));
1109566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(size, &H));
1119566063dSJacob Faibussowitsch PetscCall(PCGetOperators(pc, &A, NULL));
1129566063dSJacob Faibussowitsch PetscCall(MatGetLocalSize(A, &m, NULL));
1139566063dSJacob Faibussowitsch PetscCall(MatGetSize(A, &M, NULL));
114e53e0a0dSJakub Kruzik Mdefl = M;
115e53e0a0dSJakub Kruzik Ndefl = M;
116e53e0a0dSJakub Kruzik for (i = 0; i < size; i++) {
117e53e0a0dSJakub Kruzik if (Mdefl % 2) {
11820cd032fSJakub Kruzik if (trunc) Mdefl = (PetscInt)PetscCeilReal(Mdefl / 2.);
11920cd032fSJakub Kruzik else Mdefl = (PetscInt)PetscFloorReal((ncoeffs + Mdefl - 1) / 2.);
12020cd032fSJakub Kruzik } else Mdefl = Mdefl / 2;
1219566063dSJacob Faibussowitsch PetscCall(PCDeflationCreateSpaceWave(comm, PETSC_DECIDE, m, Mdefl, Ndefl, ncoeffs, coeffs, trunc, &H[i]));
1229566063dSJacob Faibussowitsch PetscCall(MatGetLocalSize(H[i], &m, NULL));
123e53e0a0dSJakub Kruzik Ndefl = Mdefl;
124e53e0a0dSJakub Kruzik }
1259566063dSJacob Faibussowitsch PetscCall(MatCreateComposite(comm, size, H, &defl));
1269566063dSJacob Faibussowitsch PetscCall(MatCompositeSetType(defl, MAT_COMPOSITE_MULTIPLICATIVE));
127e53e0a0dSJakub Kruzik *W = defl;
128e53e0a0dSJakub Kruzik
12948a46eb9SPierre Jolivet for (i = 0; i < size; i++) PetscCall(MatDestroy(&H[i]));
1309566063dSJacob Faibussowitsch PetscCall(PetscFree(H));
1313ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
132e53e0a0dSJakub Kruzik }
133e53e0a0dSJakub Kruzik
PCDeflationGetSpaceAggregation(PC pc,Mat * W)13466976f2fSJacob Faibussowitsch static PetscErrorCode PCDeflationGetSpaceAggregation(PC pc, Mat *W)
135d71ae5a4SJacob Faibussowitsch {
136e53e0a0dSJakub Kruzik Mat A, defl;
1377b3faf33SJakub Kruzik PetscInt i, ilo, ihi, *Iidx, M;
1387b3faf33SJakub Kruzik PetscMPIInt m;
1399e56ec8aSJakub Kruzik PetscScalar *col;
140e53e0a0dSJakub Kruzik MPI_Comm comm;
141e53e0a0dSJakub Kruzik
142e53e0a0dSJakub Kruzik PetscFunctionBegin;
1439566063dSJacob Faibussowitsch PetscCall(PCGetOperators(pc, &A, NULL));
1449566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRangeColumn(A, &ilo, &ihi));
1459566063dSJacob Faibussowitsch PetscCall(MatGetSize(A, &M, NULL));
1469566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)A, &comm));
1479566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(comm, &m));
1489566063dSJacob Faibussowitsch PetscCall(MatCreate(comm, &defl));
1499566063dSJacob Faibussowitsch PetscCall(MatSetSizes(defl, ihi - ilo, 1, M, m));
1509566063dSJacob Faibussowitsch PetscCall(MatSetUp(defl));
1519566063dSJacob Faibussowitsch PetscCall(MatSeqAIJSetPreallocation(defl, 1, NULL));
1529566063dSJacob Faibussowitsch PetscCall(MatMPIAIJSetPreallocation(defl, 1, NULL, 0, NULL));
1539566063dSJacob Faibussowitsch PetscCall(MatSetOption(defl, MAT_NEW_NONZERO_ALLOCATION_ERR, PETSC_TRUE));
1549566063dSJacob Faibussowitsch PetscCall(MatSetOption(defl, MAT_NO_OFF_PROC_ENTRIES, PETSC_TRUE));
155e53e0a0dSJakub Kruzik
1569566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(ihi - ilo, &col, ihi - ilo, &Iidx));
157e53e0a0dSJakub Kruzik for (i = ilo; i < ihi; i++) {
158e53e0a0dSJakub Kruzik Iidx[i - ilo] = i;
159e53e0a0dSJakub Kruzik col[i - ilo] = 1;
160e53e0a0dSJakub Kruzik }
1619566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(comm, &m));
1627b3faf33SJakub Kruzik i = m;
1639566063dSJacob Faibussowitsch PetscCall(MatSetValues(defl, ihi - ilo, Iidx, 1, &i, col, INSERT_VALUES));
164e53e0a0dSJakub Kruzik
1659566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(defl, MAT_FINAL_ASSEMBLY));
1669566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(defl, MAT_FINAL_ASSEMBLY));
167e53e0a0dSJakub Kruzik
1689566063dSJacob Faibussowitsch PetscCall(PetscFree2(col, Iidx));
169e53e0a0dSJakub Kruzik *W = defl;
1703ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
171e53e0a0dSJakub Kruzik }
172e53e0a0dSJakub Kruzik
PCDeflationComputeSpace(PC pc)173d71ae5a4SJacob Faibussowitsch PetscErrorCode PCDeflationComputeSpace(PC pc)
174d71ae5a4SJacob Faibussowitsch {
175e53e0a0dSJakub Kruzik Mat defl;
176e53e0a0dSJakub Kruzik PetscBool transp = PETSC_TRUE;
177e53e0a0dSJakub Kruzik PC_Deflation *def = (PC_Deflation *)pc->data;
178e53e0a0dSJakub Kruzik
179e53e0a0dSJakub Kruzik PetscFunctionBegin;
1801fdb00f9SJakub Kruzik PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
18163a3b9bcSJacob Faibussowitsch PetscCheck(def->spacesize >= 1, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONG, "Wrong PCDeflation space size specified: %" PetscInt_FMT, def->spacesize);
182e53e0a0dSJakub Kruzik switch (def->spacetype) {
183e53e0a0dSJakub Kruzik case PC_DEFLATION_SPACE_HAAR:
184e53e0a0dSJakub Kruzik transp = PETSC_FALSE;
1859371c9d4SSatish Balay PetscCall(PCDeflationGetSpaceHaar(pc, &defl, def->spacesize));
1869371c9d4SSatish Balay break;
187d71ae5a4SJacob Faibussowitsch case PC_DEFLATION_SPACE_DB2:
188d71ae5a4SJacob Faibussowitsch PetscCall(PCDeflationGetSpaceWave(pc, &defl, def->spacesize, 2, db2, PetscNot(def->extendsp)));
189d71ae5a4SJacob Faibussowitsch break;
190d71ae5a4SJacob Faibussowitsch case PC_DEFLATION_SPACE_DB4:
191d71ae5a4SJacob Faibussowitsch PetscCall(PCDeflationGetSpaceWave(pc, &defl, def->spacesize, 4, db4, PetscNot(def->extendsp)));
192d71ae5a4SJacob Faibussowitsch break;
193d71ae5a4SJacob Faibussowitsch case PC_DEFLATION_SPACE_DB8:
194d71ae5a4SJacob Faibussowitsch PetscCall(PCDeflationGetSpaceWave(pc, &defl, def->spacesize, 8, db8, PetscNot(def->extendsp)));
195d71ae5a4SJacob Faibussowitsch break;
196d71ae5a4SJacob Faibussowitsch case PC_DEFLATION_SPACE_DB16:
197d71ae5a4SJacob Faibussowitsch PetscCall(PCDeflationGetSpaceWave(pc, &defl, def->spacesize, 16, db16, PetscNot(def->extendsp)));
198d71ae5a4SJacob Faibussowitsch break;
199d71ae5a4SJacob Faibussowitsch case PC_DEFLATION_SPACE_BIORTH22:
200d71ae5a4SJacob Faibussowitsch PetscCall(PCDeflationGetSpaceWave(pc, &defl, def->spacesize, 6, biorth22, PetscNot(def->extendsp)));
201d71ae5a4SJacob Faibussowitsch break;
202d71ae5a4SJacob Faibussowitsch case PC_DEFLATION_SPACE_MEYER:
203d71ae5a4SJacob Faibussowitsch PetscCall(PCDeflationGetSpaceWave(pc, &defl, def->spacesize, 62, meyer, PetscNot(def->extendsp)));
204d71ae5a4SJacob Faibussowitsch break;
205e53e0a0dSJakub Kruzik case PC_DEFLATION_SPACE_AGGREGATION:
206e53e0a0dSJakub Kruzik transp = PETSC_FALSE;
2079371c9d4SSatish Balay PetscCall(PCDeflationGetSpaceAggregation(pc, &defl));
2089371c9d4SSatish Balay break;
209d71ae5a4SJacob Faibussowitsch default:
210d71ae5a4SJacob Faibussowitsch SETERRQ(PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONG, "Wrong PCDeflationSpaceType specified");
211e53e0a0dSJakub Kruzik }
212e53e0a0dSJakub Kruzik
2139566063dSJacob Faibussowitsch PetscCall(PCDeflationSetSpace(pc, defl, transp));
2149566063dSJacob Faibussowitsch PetscCall(MatDestroy(&defl));
2153ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
216e53e0a0dSJakub Kruzik }
217