Actual source code: pmat.c

  1: #include "tao_general.h"

  3: #ifdef TAO_USE_PETSC
  4: #include "petscmat.h"


  7: int MatD_Fischer(Mat, Vec, Vec, Vec, Vec, Vec, Vec, Vec, Vec);
  8: int MatD_SFischer(Mat, Vec, Vec, Vec, Vec, double, Vec, Vec, Vec, Vec, Vec);

 12: inline static PetscScalar Fischer(PetscScalar a, PetscScalar b)
 13: {

 15: #ifdef TODD

 17:   if (PetscAbsScalar(a) > PetscAbsScalar(b)) {
 18:     return sqrt(a*a + b*b) - a - b;
 19:   }
 20:   return sqrt(a*a + b*b) - b - a;

 22: #else

 24:    // Method suggested by Bob Vanderbei

 26:    if (a + b <= 0) {
 27:      return sqrt(a*a + b*b) - (a + b);
 28:    }
 29:    return -2.0*a*b / (sqrt(a*a + b*b) + (a + b));

 31: #endif

 33: }

 37: inline static PetscScalar Norm(PetscScalar a, PetscScalar b)
 38: {
 39:   return sqrt(a*a + b*b);
 40: }

 44: int MatD_Fischer(Mat m, Vec X, Vec F, Vec L, Vec U,
 45:                  Vec T1, Vec T2, Vec DA, Vec DB)
 46: {
 47:   PetscScalar *x, *f, *l, *u, *da, *db, *t1, *t2;
 48:   PetscScalar ai, bi, ci, di, ei;
 49:   int info, i, n;
 50:   int low[8], high[8];


 62:   info = VecGetOwnershipRange(X, low, high); CHKERRQ(info);
 63:   info = VecGetOwnershipRange(F, low + 1, high + 1); CHKERRQ(info);
 64:   info = VecGetOwnershipRange(L, low + 2, high + 2); CHKERRQ(info);
 65:   info = VecGetOwnershipRange(U, low + 3, high + 3); CHKERRQ(info);
 66:   info = VecGetOwnershipRange(DA, low + 4, high + 4); CHKERRQ(info);
 67:   info = VecGetOwnershipRange(DB, low + 5, high + 5); CHKERRQ(info);
 68:   info = VecGetOwnershipRange(T1, low + 6, high + 6); CHKERRQ(info);
 69:   info = VecGetOwnershipRange(T2, low + 7, high + 7); CHKERRQ(info);

 71:   for (i = 1; i < 8; i++) {
 72:     if (low[0] != low[i] || high[0] != high[i])
 73:       SETERRQ(1,"Vectors must be identically loaded over processors");
 74:   }

 76:   info = VecGetArray(X, &x); CHKERRQ(info);
 77:   info = VecGetArray(F, &f); CHKERRQ(info);
 78:   info = VecGetArray(L, &l); CHKERRQ(info);
 79:   info = VecGetArray(U, &u); CHKERRQ(info);
 80:   info = VecGetArray(DA, &da); CHKERRQ(info);
 81:   info = VecGetArray(DB, &db); CHKERRQ(info);
 82:   info = VecGetArray(T1, &t1); CHKERRQ(info);

 84:   info = VecGetLocalSize(X, &n); CHKERRQ(info);

 86:   for (i = 0; i < n; i++) {
 87:     da[i] = 0;
 88:     db[i] = 0;
 89:     t1[i] = 0;

 91:     if (PetscAbsScalar(f[i]) <= TAO_EPSILON) {
 92:       if (l[i] > -TAO_INFINITY && PetscAbsScalar(x[i] - l[i]) <= TAO_EPSILON) {
 93:         t1[i] = 1;
 94:         da[i] = 1;
 95:       }

 97:       if (u[i] <  TAO_INFINITY && PetscAbsScalar(u[i] - x[i]) <= TAO_EPSILON) {
 98:         t1[i] = 1;
 99:         db[i] = 1;
100:       }
101:     }
102:   }

104:   info = VecRestoreArray(T1, &t1); CHKERRQ(info);
105:   info = MatMult(m, T1, T2); CHKERRQ(info);
106:   info = VecGetArray(T2, &t2); CHKERRQ(info);
107:   
108:   for (i = 0; i < n; i++) {
109:     if ((l[i] <= -TAO_INFINITY) && (u[i] >= TAO_INFINITY)) {
110:       da[i] = 0;
111:       db[i] = -1;
112:     } 
113:     else if (l[i] <= -TAO_INFINITY) {
114:       if (db[i] >= 1) {
115:         ai = Norm(1, t2[i]);

117:         da[i] = -1/ai - 1;
118:         db[i] = -t2[i]/ai - 1;
119:       } 
120:       else {
121:         bi = u[i] - x[i];
122:         ai = Norm(bi, f[i]);
123:         ai = PetscMax(TAO_EPSILON, ai);

125:         da[i] = bi / ai - 1;
126:         db[i] = -f[i] / ai - 1;
127:       }
128:     } 
129:     else if (u[i] >=  TAO_INFINITY) {
130:       if (da[i] >= 1) {
131:         ai = Norm(1, t2[i]);

133:         da[i] = 1 / ai - 1;
134:         db[i] = t2[i] / ai - 1;
135:       } 
136:       else {
137:         bi = x[i] - l[i];
138:         ai = Norm(bi, f[i]);
139:         ai = PetscMax(TAO_EPSILON, ai);

141:         da[i] = bi / ai - 1;
142:         db[i] = f[i] / ai - 1;
143:       }
144:     } 
145:     else if (l[i] == u[i]) {
146:       da[i] = -1;
147:       db[i] = 0;
148:     } 
149:     else {
150:       if (db[i] >= 1) {
151:         ai = Norm(1, t2[i]);

153:         ci = 1 / ai + 1;
154:         di = t2[i] / ai + 1;
155:       } 
156:       else {
157:         bi = x[i] - u[i];
158:         ai = Norm(bi, f[i]);
159:         ai = PetscMax(TAO_EPSILON, ai);

161:         ci = bi / ai + 1;
162:         di = f[i] / ai + 1;
163:       }

165:       if (da[i] >= 1) {
166:         bi = ci + di*t2[i];
167:         ai = Norm(1, bi);
168:         
169:         bi = bi / ai - 1;
170:         ai = 1 / ai - 1;
171:       } 
172:       else {
173:         ei = Fischer(u[i] - x[i], -f[i]);
174:         ai = Norm(x[i] - l[i], ei);
175:         ai = PetscMax(TAO_EPSILON, ai);

177:         bi = ei / ai - 1;
178:         ai = (x[i] - l[i]) / ai - 1;
179:       }

181:       da[i] = ai + bi*ci;
182:       db[i] = bi*di;
183:     }
184:   }

186:   info = VecRestoreArray(DA, &da); CHKERRQ(info);
187:   info = VecRestoreArray(DB, &db); CHKERRQ(info);

189:   info = VecRestoreArray(X, &x); CHKERRQ(info);
190:   info = VecRestoreArray(F, &f); CHKERRQ(info);
191:   info = VecRestoreArray(L, &l); CHKERRQ(info);
192:   info = VecRestoreArray(U, &u); CHKERRQ(info);
193:   info = VecRestoreArray(T2, &t2); CHKERRQ(info);
194:   return(0);
195: }

199: inline static PetscScalar SFischer(PetscScalar a, PetscScalar b, PetscScalar c)
200: {

202: #ifdef TODD

204:   if (PetscAbsScalar(a) > PetscAbsScalar(b)) {
205:     return sqrt(a*a + b*b + 2.0*c*c) - a - b;
206:   }
207:   return sqrt(a*a + b*b + 2.0*c*c) - b - a;

209: #else

211:    // Method suggested by Bob Vanderbei

213:    if (a + b <= 0) {
214:      return sqrt(a*a + b*b + 2.0*c*c) - (a + b);
215:    }
216:    return 2.0*(c*c - a*b) / (sqrt(a*a + b*b + 2.0*c*c) + (a + b));

218: #endif

220: }

224: inline static PetscScalar SNorm(PetscScalar a, PetscScalar b, PetscScalar c)
225: {
226:   return sqrt(a*a + b*b + 2.0*c*c);
227: }

231: int MatD_SFischer(Mat m, Vec X, Vec F, Vec L, Vec U, double mu,
232:                   Vec T1, Vec T2, Vec DA, Vec DB, Vec DM)
233: {
234:   PetscScalar *x, *f, *l, *u, *da, *db, *dm;
235:   PetscScalar ai, bi, ci, di, ei, fi;
236:   int info, i, n;
237:   int low[7], high[7];


248:   info = VecGetOwnershipRange(X, low, high); CHKERRQ(info);
249:   info = VecGetOwnershipRange(F, low + 1, high + 1); CHKERRQ(info);
250:   info = VecGetOwnershipRange(L, low + 2, high + 2); CHKERRQ(info);
251:   info = VecGetOwnershipRange(U, low + 3, high + 3); CHKERRQ(info);
252:   info = VecGetOwnershipRange(DA, low + 4, high + 4); CHKERRQ(info);
253:   info = VecGetOwnershipRange(DB, low + 5, high + 5); CHKERRQ(info);
254:   info = VecGetOwnershipRange(DM, low + 6, high + 6); CHKERRQ(info);

256:   for (i = 1; i < 7; i++) {
257:     if (low[0] != low[i] || high[0] != high[i])
258:       SETERRQ(1,"Vectors must be identically loaded over processors");
259:   }

261:   info = VecGetArray(X, &x); CHKERRQ(info);
262:   info = VecGetArray(F, &f); CHKERRQ(info);
263:   info = VecGetArray(L, &l); CHKERRQ(info);
264:   info = VecGetArray(U, &u); CHKERRQ(info);
265:   info = VecGetArray(DA, &da); CHKERRQ(info);
266:   info = VecGetArray(DB, &db); CHKERRQ(info);
267:   info = VecGetArray(DM, &dm); CHKERRQ(info);

269:   info = VecGetLocalSize(X, &n); CHKERRQ(info);

271:   for (i = 0; i < n; i++) {
272:     if ((l[i] <= -TAO_INFINITY) && (u[i] >= TAO_INFINITY)) {
273:       da[i] = -mu;
274:       db[i] = -1;
275:       dm[i] = -x[i];
276:     } 
277:     else if (l[i] <= -TAO_INFINITY) {
278:       bi = u[i] - x[i];
279:       ai = SNorm(bi, f[i], mu);
280:       ai = PetscMax(TAO_EPSILON, ai);

282:       da[i] = bi / ai - 1;
283:       db[i] = -f[i] / ai - 1;
284:       dm[i] = 2.0 * mu / ai;
285:     } 
286:     else if (u[i] >=  TAO_INFINITY) {
287:       bi = x[i] - l[i];
288:       ai = SNorm(bi, f[i], mu);
289:       ai = PetscMax(TAO_EPSILON, ai);

291:       da[i] = bi / ai - 1;
292:       db[i] = f[i] / ai - 1;
293:       dm[i] = 2.0 * mu / ai;
294:     } 
295:     else if (l[i] == u[i]) {
296:       da[i] = -1;
297:       db[i] = 0;
298:       dm[i] = 0;
299:     } 
300:     else {
301:       bi = x[i] - u[i];
302:       ai = SNorm(bi, f[i], mu);
303:       ai = PetscMax(TAO_EPSILON, ai);

305:       ci = bi / ai + 1;
306:       di = f[i] / ai + 1;
307:       fi = 2.0 * mu / ai;

309:       ei = SFischer(u[i] - x[i], -f[i], mu);
310:       ai = SNorm(x[i] - l[i], ei, mu);
311:       ai = PetscMax(TAO_EPSILON, ai);

313:       bi = ei / ai - 1;
314:       ai = (x[i] - l[i]) / ai - 1;

316:       da[i] = ai + bi*ci;
317:       db[i] = bi*di;
318:       dm[i] = ei + bi*fi;
319:     }
320:   }

322:   info = VecRestoreArray(X, &x); CHKERRQ(info);
323:   info = VecRestoreArray(F, &f); CHKERRQ(info);
324:   info = VecRestoreArray(L, &l); CHKERRQ(info);
325:   info = VecRestoreArray(U, &u); CHKERRQ(info);
326:   info = VecRestoreArray(DA, &da); CHKERRQ(info);
327:   info = VecRestoreArray(DB, &db); CHKERRQ(info);
328:   info = VecRestoreArray(DM, &dm); CHKERRQ(info);
329:   return(0);
330: }

332: #endif