Actual source code: agmresorthog.c
1: #define PETSCKSP_DLL
3: #include <../src/ksp/ksp/impls/gmres/agmres/agmresimpl.h>
4: /*
5: * This file implements the RODDEC algorithm : its purpose is to orthogonalize a set of vectors distributed across several processes. These processes are organized in a virtual ring.
6: * References : [1] Sidje, Roger B. Alternatives for parallel Krylov subspace basis computation. Numer. Linear Algebra Appl. 4 (1997), no. 4, 305-331
7: *
8: *
9: * initial author R. B. SIDJE,
10: * modified : G.-A Atenekeng-Kahou, 2008
11: * modified : D. NUENTSA WAKAM, 2011
12: *
13: */
15: /*
16: * Take the processes that own the vectors and Organize them on a virtual ring.
17: */
18: static PetscErrorCode KSPAGMRESRoddecGivens(PetscReal *, PetscReal *, PetscReal *, PetscInt);
20: PetscErrorCode KSPAGMRESRoddecInitNeighboor(KSP ksp)
21: {
22: MPI_Comm comm;
23: KSP_AGMRES *agmres = (KSP_AGMRES *)(ksp->data);
24: PetscMPIInt First, Last, rank, size;
26: PetscObjectGetComm((PetscObject)agmres->vecs[0], &comm);
27: MPI_Comm_rank(comm, &rank);
28: MPI_Comm_size(comm, &size);
29: MPIU_Allreduce(&rank, &First, 1, MPI_INT, MPI_MIN, comm);
30: MPIU_Allreduce(&rank, &Last, 1, MPI_INT, MPI_MAX, comm);
32: if ((rank != Last) && (rank == 0)) {
33: agmres->Ileft = rank - 1;
34: agmres->Iright = rank + 1;
35: } else {
36: if (rank == Last) {
37: agmres->Ileft = rank - 1;
38: agmres->Iright = First;
39: } else {
40: {
41: agmres->Ileft = Last;
42: agmres->Iright = rank + 1;
43: }
44: }
45: }
46: agmres->rank = rank;
47: agmres->size = size;
48: agmres->First = First;
49: agmres->Last = Last;
50: return 0;
51: }
53: static PetscErrorCode KSPAGMRESRoddecGivens(PetscReal *c, PetscReal *s, PetscReal *r, PetscInt make_r)
54: {
55: PetscReal a, b, t;
57: if (make_r == 1) {
58: a = *c;
59: b = *s;
60: if (b == 0.e0) {
61: *c = 1.e0;
62: *s = 0.e0;
63: } else {
64: if (PetscAbsReal(b) > PetscAbsReal(a)) {
65: t = -a / b;
66: *s = 1.e0 / PetscSqrtReal(1.e0 + t * t);
67: *c = (*s) * t;
68: } else {
69: t = -b / a;
70: *c = 1.e0 / PetscSqrtReal(1.e0 + t * t);
71: *s = (*c) * t;
72: }
73: }
74: if (*c == 0.e0) {
75: *r = 1.e0;
76: } else {
77: if (PetscAbsReal(*s) < PetscAbsReal(*c)) {
78: *r = PetscSign(*c) * (*s) / 2.e0;
79: } else {
80: *r = PetscSign(*s) * 2.e0 / (*c);
81: }
82: }
83: }
85: if (*r == 1.e0) {
86: *c = 0.e0;
87: *s = 1.e0;
88: } else {
89: if (PetscAbsReal(*r) < 1.e0) {
90: *s = 2.e0 * (*r);
91: *c = PetscSqrtReal(1.e0 - (*s) * (*s));
92: } else {
93: *c = 2.e0 / (*r);
94: *s = PetscSqrtReal(1.e0 - (*c) * (*c));
95: }
96: }
97: return 0;
98: }
100: /*
101: * Compute the QR factorization of the Krylov basis vectors
102: * Input :
103: * - the vectors are available in agmres->vecs (alias VEC_V)
104: * - nvec : the number of vectors
105: * Output :
106: * - agmres->Qloc : product of elementary reflectors for the QR of the (local part) of the vectors.
107: * - agmres->sgn : Sign of the rotations
108: * - agmres->tloc : scalar factors of the elementary reflectors.
110: */
111: PetscErrorCode KSPAGMRESRoddec(KSP ksp, PetscInt nvec)
112: {
113: KSP_AGMRES *agmres = (KSP_AGMRES *)ksp->data;
114: MPI_Comm comm;
115: PetscScalar *Qloc = agmres->Qloc;
116: PetscScalar *sgn = agmres->sgn;
117: PetscScalar *tloc = agmres->tloc;
118: PetscReal *wbufptr = agmres->wbufptr;
119: PetscMPIInt rank = agmres->rank;
120: PetscMPIInt First = agmres->First;
121: PetscMPIInt Last = agmres->Last;
122: PetscBLASInt pas, len, bnloc, bpos;
123: PetscInt nloc, d, i, j, k;
124: PetscInt pos;
125: PetscReal c, s, rho, Ajj, val, tt, old;
126: PetscScalar *col;
127: MPI_Status status;
128: PetscBLASInt N = MAXKSPSIZE + 1;
130: PetscObjectGetComm((PetscObject)ksp, &comm);
131: PetscLogEventBegin(KSP_AGMRESRoddec, ksp, 0, 0, 0);
132: PetscArrayzero(agmres->Rloc, N * N);
133: /* check input arguments */
135: VecGetLocalSize(VEC_V(0), &nloc);
136: PetscBLASIntCast(nloc, &bnloc);
138: pas = 1;
139: /* Copy the vectors of the basis */
140: for (j = 0; j < nvec; j++) {
141: VecGetArray(VEC_V(j), &col);
142: PetscCallBLAS("BLAScopy", BLAScopy_(&bnloc, col, &pas, &Qloc[j * nloc], &pas));
143: VecRestoreArray(VEC_V(j), &col);
144: }
145: /* Each process performs a local QR on its own block */
146: for (j = 0; j < nvec; j++) {
147: len = nloc - j;
148: Ajj = Qloc[j * nloc + j];
149: PetscCallBLAS("BLASnrm2", rho = -PetscSign(Ajj) * BLASnrm2_(&len, &(Qloc[j * nloc + j]), &pas));
150: if (rho == 0.0) tloc[j] = 0.0;
151: else {
152: tloc[j] = (Ajj - rho) / rho;
153: len = len - 1;
154: val = 1.0 / (Ajj - rho);
155: PetscCallBLAS("BLASscal", BLASscal_(&len, &val, &(Qloc[j * nloc + j + 1]), &pas));
156: Qloc[j * nloc + j] = 1.0;
157: len = len + 1;
158: for (k = j + 1; k < nvec; k++) {
159: PetscCallBLAS("BLASdot", tt = tloc[j] * BLASdot_(&len, &(Qloc[j * nloc + j]), &pas, &(Qloc[k * nloc + j]), &pas));
160: PetscCallBLAS("BLASaxpy", BLASaxpy_(&len, &tt, &(Qloc[j * nloc + j]), &pas, &(Qloc[k * nloc + j]), &pas));
161: }
162: Qloc[j * nloc + j] = rho;
163: }
164: }
165: /* annihilate undesirable Rloc, diagonal by diagonal*/
166: for (d = 0; d < nvec; d++) {
167: len = nvec - d;
168: if (rank == First) {
169: PetscCallBLAS("BLAScopy", BLAScopy_(&len, &(Qloc[d * nloc + d]), &bnloc, &(wbufptr[d]), &pas));
170: MPI_Send(&(wbufptr[d]), len, MPIU_SCALAR, rank + 1, agmres->tag, comm);
171: } else {
172: MPI_Recv(&(wbufptr[d]), len, MPIU_SCALAR, rank - 1, agmres->tag, comm, &status);
173: /* Elimination of Rloc(1,d)*/
174: c = wbufptr[d];
175: s = Qloc[d * nloc];
176: KSPAGMRESRoddecGivens(&c, &s, &rho, 1);
177: /* Apply Givens Rotation*/
178: for (k = d; k < nvec; k++) {
179: old = wbufptr[k];
180: wbufptr[k] = c * old - s * Qloc[k * nloc];
181: Qloc[k * nloc] = s * old + c * Qloc[k * nloc];
182: }
183: Qloc[d * nloc] = rho;
184: if (rank != Last) MPI_Send(&(wbufptr[d]), len, MPIU_SCALAR, rank + 1, agmres->tag, comm);
185: /* zero-out the d-th diagonal of Rloc ...*/
186: for (j = d + 1; j < nvec; j++) {
187: /* elimination of Rloc[i][j]*/
188: i = j - d;
189: c = Qloc[j * nloc + i - 1];
190: s = Qloc[j * nloc + i];
191: KSPAGMRESRoddecGivens(&c, &s, &rho, 1);
192: for (k = j; k < nvec; k++) {
193: old = Qloc[k * nloc + i - 1];
194: Qloc[k * nloc + i - 1] = c * old - s * Qloc[k * nloc + i];
195: Qloc[k * nloc + i] = s * old + c * Qloc[k * nloc + i];
196: }
197: Qloc[j * nloc + i] = rho;
198: }
199: if (rank == Last) {
200: PetscCallBLAS("BLAScopy", BLAScopy_(&len, &(wbufptr[d]), &pas, RLOC(d, d), &N));
201: for (k = d + 1; k < nvec; k++) *RLOC(k, d) = 0.0;
202: }
203: }
204: }
206: if (rank == Last) {
207: for (d = 0; d < nvec; d++) {
208: pos = nvec - d;
209: PetscBLASIntCast(pos, &bpos);
210: sgn[d] = PetscSign(*RLOC(d, d));
211: PetscCallBLAS("BLASscal", BLASscal_(&bpos, &(sgn[d]), RLOC(d, d), &N));
212: }
213: }
214: /* BroadCast Rloc to all other processes
215: * NWD : should not be needed
216: */
217: MPI_Bcast(agmres->Rloc, N * N, MPIU_SCALAR, Last, comm);
218: PetscLogEventEnd(KSP_AGMRESRoddec, ksp, 0, 0, 0);
219: return 0;
220: }
222: /*
223: * Computes Out <-- Q * In where Q is the orthogonal matrix from AGMRESRoddec
224: * Input
225: * - Qloc, sgn, tloc, nvec (see AGMRESRoddec above)
226: * - In : input vector (size nvec)
227: * Output :
228: * - Out : Petsc vector (distributed as the basis vectors)
229: */
230: PetscErrorCode KSPAGMRESRodvec(KSP ksp, PetscInt nvec, PetscScalar *In, Vec Out)
231: {
232: KSP_AGMRES *agmres = (KSP_AGMRES *)ksp->data;
233: MPI_Comm comm;
234: PetscScalar *Qloc = agmres->Qloc;
235: PetscScalar *sgn = agmres->sgn;
236: PetscScalar *tloc = agmres->tloc;
237: PetscMPIInt rank = agmres->rank;
238: PetscMPIInt First = agmres->First, Last = agmres->Last;
239: PetscMPIInt Iright = agmres->Iright, Ileft = agmres->Ileft;
240: PetscScalar *y, *zloc;
241: PetscInt nloc, d, len, i, j;
242: PetscBLASInt bnvec, pas, blen;
243: PetscInt dpt;
244: PetscReal c, s, rho, zp, zq, yd = 0.0, tt;
245: MPI_Status status;
247: PetscBLASIntCast(nvec, &bnvec);
248: PetscObjectGetComm((PetscObject)ksp, &comm);
249: pas = 1;
250: VecGetLocalSize(VEC_V(0), &nloc);
251: PetscMalloc1(nvec, &y);
252: PetscArraycpy(y, In, nvec);
253: VecGetArray(Out, &zloc);
255: if (rank == Last) {
256: for (i = 0; i < nvec; i++) y[i] = sgn[i] * y[i];
257: }
258: for (i = 0; i < nloc; i++) zloc[i] = 0.0;
259: if (agmres->size == 1) PetscCallBLAS("BLAScopy", BLAScopy_(&bnvec, y, &pas, &(zloc[0]), &pas));
260: else {
261: for (d = nvec - 1; d >= 0; d--) {
262: if (rank == First) {
263: MPI_Recv(&(zloc[d]), 1, MPIU_SCALAR, Iright, agmres->tag, comm, &status);
264: } else {
265: for (j = nvec - 1; j >= d + 1; j--) {
266: i = j - d;
267: KSPAGMRESRoddecGivens(&c, &s, &(Qloc[j * nloc + i]), 0);
268: zp = zloc[i - 1];
269: zq = zloc[i];
270: zloc[i - 1] = c * zp + s * zq;
271: zloc[i] = -s * zp + c * zq;
272: }
273: KSPAGMRESRoddecGivens(&c, &s, &(Qloc[d * nloc]), 0);
274: if (rank == Last) {
275: zp = y[d];
276: zq = zloc[0];
277: y[d] = c * zp + s * zq;
278: zloc[0] = -s * zp + c * zq;
279: MPI_Send(&(y[d]), 1, MPIU_SCALAR, Ileft, agmres->tag, comm);
280: } else {
281: MPI_Recv(&yd, 1, MPIU_SCALAR, Iright, agmres->tag, comm, &status);
282: zp = yd;
283: zq = zloc[0];
284: yd = c * zp + s * zq;
285: zloc[0] = -s * zp + c * zq;
286: MPI_Send(&yd, 1, MPIU_SCALAR, Ileft, agmres->tag, comm);
287: }
288: }
289: }
290: }
291: for (j = nvec - 1; j >= 0; j--) {
292: dpt = j * nloc + j;
293: if (tloc[j] != 0.0) {
294: len = nloc - j;
295: PetscBLASIntCast(len, &blen);
296: rho = Qloc[dpt];
297: Qloc[dpt] = 1.0;
298: tt = tloc[j] * (BLASdot_(&blen, &(Qloc[dpt]), &pas, &(zloc[j]), &pas));
299: PetscCallBLAS("BLASaxpy", BLASaxpy_(&blen, &tt, &(Qloc[dpt]), &pas, &(zloc[j]), &pas));
300: Qloc[dpt] = rho;
301: }
302: }
303: VecRestoreArray(Out, &zloc);
304: PetscFree(y);
305: return 0;
306: }