Actual source code: pvec.c
1: #include "tao_general.h"
2: #include "private/vecimpl.h" /*I "tao_solver.h" I*/
3: #include "petscksp.h"
9: int VecCompare(Vec V1,Vec V2, PetscTruth *flg){
10: int info;
11: int n1,n2,N1,N2;
16: info = VecGetSize(V1,&N1);CHKERRQ(info);
17: info = VecGetSize(V2,&N2);CHKERRQ(info);
18: info = VecGetLocalSize(V1,&n1);CHKERRQ(info);
19: info = VecGetLocalSize(V2,&n2);CHKERRQ(info);
20: if (N1==N2 && n1==n2)
21: *flg=PETSC_TRUE;
22: else
23: *flg=PETSC_FALSE;
25: return(0);
26: }
28: /* ---------------------------------------------------------- */
31: int VecMedian(Vec Vec1, Vec Vec2, Vec Vec3, Vec VMedian)
32: {
33: int i,n,low1,low2,low3,low4,high1,high2,high3,high4,info;
34: PetscScalar *v1,*v2,*v3,*vmed;
42: if (Vec1==Vec2 || Vec1==Vec3){
43: info=VecCopy(Vec1,VMedian); CHKERRQ(info);
44: return(0);
45: }
46: if (Vec2==Vec3){
47: info=VecCopy(Vec2,VMedian); CHKERRQ(info);
48: return(0);
49: }
57: info = VecGetOwnershipRange(Vec1, &low1, &high1); CHKERRQ(info);
58: info = VecGetOwnershipRange(Vec2, &low2, &high2); CHKERRQ(info);
59: info = VecGetOwnershipRange(Vec3, &low3, &high3); CHKERRQ(info);
60: info = VecGetOwnershipRange(VMedian, &low4, &high4); CHKERRQ(info);
61: if ( low1!= low2 || low1!= low3 || low1!= low4 ||
62: high1!= high2 || high1!= high3 || high1!= high4)
63: SETERRQ(1,"InCompatible vector local lengths");
65: info = VecGetArray(Vec1,&v1); CHKERRQ(info);
66: info = VecGetArray(Vec2,&v2); CHKERRQ(info);
67: info = VecGetArray(Vec3,&v3); CHKERRQ(info);
69: if ( VMedian != Vec1 && VMedian != Vec2 && VMedian != Vec3){
70: info = VecGetArray(VMedian,&vmed); CHKERRQ(info);
71: } else if ( VMedian==Vec1 ){
72: vmed=v1;
73: } else if ( VMedian==Vec2 ){
74: vmed=v2;
75: } else {
76: vmed=v3;
77: }
79: info=VecGetLocalSize(Vec1,&n); CHKERRQ(info);
81: for (i=0;i<n;i++){
82: vmed[i]=PetscMax(PetscMax(PetscMin(v1[i],v2[i]),PetscMin(v1[i],v3[i])),PetscMin(v2[i],v3[i]));
83: }
85: info = VecRestoreArray(Vec1,&v1); CHKERRQ(info);
86: info = VecRestoreArray(Vec2,&v2); CHKERRQ(info);
87: info = VecRestoreArray(Vec3,&v2); CHKERRQ(info);
88:
89: if (VMedian!=Vec1 && VMedian != Vec2 && VMedian != Vec3){
90: info = VecRestoreArray(VMedian,&vmed); CHKERRQ(info);
91: }
93: return(0);
94: }
96: /* ----------------------------------------------------------
99: int VecPointwiseMin(Vec Vec1, Vec Vec2, Vec VMin)
100: {
101: int i,n,low1,low2,low3,high1,high2,high3,info;
102: PetscScalar *v1,*v2,*vmin;
109: if (Vec1==Vec2){
110: info=VecCopy(Vec1,VMin); CHKERRQ(info);
111: return(0);
112: }
117: info = VecGetOwnershipRange(Vec1, &low1, &high1); CHKERRQ(info);
118: info = VecGetOwnershipRange(Vec2, &low2, &high2); CHKERRQ(info);
119: info = VecGetOwnershipRange(VMin, &low3, &high3); CHKERRQ(info);
121: if ( low1!= low2 || low1!= low3 || high1!= high2 || high1!= high3)
122: SETERRQ(1,"InCompatible vector local lengths");
124: info = VecGetArray(Vec1,&v1); CHKERRQ(info);
125: info = VecGetArray(Vec2,&v2); CHKERRQ(info);
127: if ( VMin != Vec1 && VMin != Vec2){
128: info = VecGetArray(VMin,&vmin); CHKERRQ(info);
129: } else if ( VMin==Vec1 ){
130: vmin=v1;
131: } else {
132: vmin =v2;
133: }
135: info=VecGetLocalSize(Vec1,&n); CHKERRQ(info);
137: for (i=0; i<n; i++){
138: vmin[i] = PetscMin(v1[i],v2[i]);
139: }
141: info = VecRestoreArray(Vec1,&v1); CHKERRQ(info);
142: info = VecRestoreArray(Vec2,&v2); CHKERRQ(info);
144: if (VMin!=Vec1 && VMin != Vec2){
145: info = VecRestoreArray(VMin,&vmin); CHKERRQ(info);
146: }
148: return(0);
149: }
150: */
151: /* ----------------------------------------------------------
154: int VecPointwiseMax(Vec Vec1, Vec Vec2, Vec VMax)
155: {
156: int i,n,low1,low2,low3,high1,high2,high3,info;
157: PetscScalar *v1,*v2,*vmax;
161: if (Vec1==Vec2){
162: info=VecCopy(Vec1,VMax); CHKERRQ(info);
163: return(0);
164: }
173: info = VecGetOwnershipRange(Vec1, &low1, &high1); CHKERRQ(info);
174: info = VecGetOwnershipRange(Vec2, &low2, &high2); CHKERRQ(info);
175: info = VecGetOwnershipRange(VMax, &low3, &high3); CHKERRQ(info);
177: if ( low1!= low2 || low1!= low3 || high1!= high2 || high1!= high3)
178: SETERRQ(1,"Vectors must be identically loaded over processors");
180: info = VecGetArray(Vec1,&v1); CHKERRQ(info);
181: info = VecGetArray(Vec2,&v2); CHKERRQ(info);
183: if ( VMax != Vec1 && VMax != Vec2){
184: info = VecGetArray(VMax,&vmax); CHKERRQ(info);
185: } else if ( VMax==Vec1 ){ vmax=v1;
186: } else { vmax =v2;
187: }
189: info=VecGetLocalSize(Vec1,&n); CHKERRQ(info);
191: for (i=0; i<n; i++){
192: vmax[i] = PetscMax(v1[i],v2[i]);
193: }
195: info = VecRestoreArray(Vec1,&v1); CHKERRQ(info);
196: info = VecRestoreArray(Vec2,&v2); CHKERRQ(info);
198: if ( VMax != Vec1 && VMax != Vec2){
199: info = VecRestoreArray(VMax,&vmax); CHKERRQ(info);
200: }
202: return(0);
203: }*/
207: inline static PetscScalar Fischer(PetscScalar a, PetscScalar b)
208: {
209: // Method suggested by Bob Vanderbei
210: if (a + b <= 0) {
211: return sqrt(a*a + b*b) - (a + b);
212: }
213: return -2.0*a*b / (sqrt(a*a + b*b) + (a + b));
214: }
218: int VecFischer(Vec X, Vec F, Vec L, Vec U, Vec FF)
219: {
220: PetscScalar *x, *f, *l, *u, *ff;
221: PetscScalar xval, fval, lval, uval;
222: int info, i, n;
223: int low[5], high[5];
232: info = VecGetOwnershipRange(X, low, high); CHKERRQ(info);
233: info = VecGetOwnershipRange(F, low + 1, high + 1); CHKERRQ(info);
234: info = VecGetOwnershipRange(L, low + 2, high + 2); CHKERRQ(info);
235: info = VecGetOwnershipRange(U, low + 3, high + 3); CHKERRQ(info);
236: info = VecGetOwnershipRange(FF, low + 4, high + 4); CHKERRQ(info);
238: for (i = 1; i < 4; ++i) {
239: if (low[0] != low[i] || high[0] != high[i])
240: SETERRQ(1,"Vectors must be identically loaded over processors");
241: }
243: info = VecGetArray(X, &x); CHKERRQ(info);
244: info = VecGetArray(F, &f); CHKERRQ(info);
245: info = VecGetArray(L, &l); CHKERRQ(info);
246: info = VecGetArray(U, &u); CHKERRQ(info);
247: info = VecGetArray(FF, &ff); CHKERRQ(info);
249: info = VecGetLocalSize(X, &n); CHKERRQ(info);
251: for (i = 0; i < n; ++i) {
252: xval = x[i]; fval = f[i];
253: lval = l[i]; uval = u[i];
255: if ((lval <= -TAO_INFINITY) && (uval >= TAO_INFINITY)) {
256: ff[i] = -fval;
257: }
258: else if (lval <= -TAO_INFINITY) {
259: ff[i] = -Fischer(uval - xval, -fval);
260: }
261: else if (uval >= TAO_INFINITY) {
262: ff[i] = Fischer(xval - lval, fval);
263: }
264: else if (lval == uval) {
265: ff[i] = lval - xval;
266: }
267: else {
268: fval = Fischer(uval - xval, -fval);
269: ff[i] = Fischer(xval - lval, fval);
270: }
271: }
272:
273: info = VecRestoreArray(X, &x); CHKERRQ(info);
274: info = VecRestoreArray(F, &f); CHKERRQ(info);
275: info = VecRestoreArray(L, &l); CHKERRQ(info);
276: info = VecRestoreArray(U, &u); CHKERRQ(info);
277: info = VecRestoreArray(FF, &ff); CHKERRQ(info);
279: return(0);
280: }
284: inline static PetscScalar SFischer(PetscScalar a, PetscScalar b, PetscScalar c)
285: {
286: // Method suggested by Bob Vanderbei
287: if (a + b <= 0) {
288: return sqrt(a*a + b*b + 2.0*c*c) - (a + b);
289: }
290: return 2.0*(c*c - a*b) / (sqrt(a*a + b*b + 2.0*c*c) + (a + b));
291: }
295: int VecSFischer(Vec X, Vec F, Vec L, Vec U, double mu, Vec FF)
296: {
297: PetscScalar *x, *f, *l, *u, *ff;
298: PetscScalar xval, fval, lval, uval;
299: int info, i, n;
300: int low[5], high[5];
309: info = VecGetOwnershipRange(X, low, high); CHKERRQ(info);
310: info = VecGetOwnershipRange(F, low + 1, high + 1); CHKERRQ(info);
311: info = VecGetOwnershipRange(L, low + 2, high + 2); CHKERRQ(info);
312: info = VecGetOwnershipRange(U, low + 3, high + 3); CHKERRQ(info);
313: info = VecGetOwnershipRange(FF, low + 4, high + 4); CHKERRQ(info);
315: for (i = 1; i < 4; ++i) {
316: if (low[0] != low[i] || high[0] != high[i])
317: SETERRQ(1,"Vectors must be identically loaded over processors");
318: }
320: info = VecGetArray(X, &x); CHKERRQ(info);
321: info = VecGetArray(F, &f); CHKERRQ(info);
322: info = VecGetArray(L, &l); CHKERRQ(info);
323: info = VecGetArray(U, &u); CHKERRQ(info);
324: info = VecGetArray(FF, &ff); CHKERRQ(info);
326: info = VecGetLocalSize(X, &n); CHKERRQ(info);
328: for (i = 0; i < n; ++i) {
329: xval = (*x++); fval = (*f++);
330: lval = (*l++); uval = (*u++);
332: if ((lval <= -TAO_INFINITY) && (uval >= TAO_INFINITY)) {
333: (*ff++) = -fval - mu*xval;
334: }
335: else if (lval <= -TAO_INFINITY) {
336: (*ff++) = -SFischer(uval - xval, -fval, mu);
337: }
338: else if (uval >= TAO_INFINITY) {
339: (*ff++) = SFischer(xval - lval, fval, mu);
340: }
341: else if (lval == uval) {
342: (*ff++) = lval - xval;
343: }
344: else {
345: fval = SFischer(uval - xval, -fval, mu);
346: (*ff++) = SFischer(xval - lval, fval, mu);
347: }
348: }
349: x -= n; f -= n; l -=n; u -= n; ff -= n;
351: info = VecRestoreArray(X, &x); CHKERRQ(info);
352: info = VecRestoreArray(F, &f); CHKERRQ(info);
353: info = VecRestoreArray(L, &l); CHKERRQ(info);
354: info = VecRestoreArray(U, &u); CHKERRQ(info);
355: info = VecRestoreArray(FF, &ff); CHKERRQ(info);
356: return(0);
357: }
361: int VecBoundProjection(Vec G, Vec X, Vec XL, Vec XU, Vec GP){
363: int i,n,info;
364: PetscScalar *xptr,*xlptr,*xuptr,*gptr,*gpptr;
365: PetscScalar xval,gpval;
367: /* Project variables at the lower and upper bound */
370: info = VecGetLocalSize(X,&n); CHKERRQ(info);
372: info=VecGetArray(X,&xptr); CHKERRQ(info);
373: info=VecGetArray(XL,&xlptr); CHKERRQ(info);
374: info=VecGetArray(XU,&xuptr); CHKERRQ(info);
375: info=VecGetArray(G,&gptr); CHKERRQ(info);
376: if (G!=GP){
377: info=VecGetArray(GP,&gpptr); CHKERRQ(info);
378: } else { gpptr=gptr; }
380: for (i=0; i<n; ++i){
381: gpval = gptr[i]; xval = xptr[i];
383: if (gpval>0 && xval<=xlptr[i]){
384: gpval = 0;
385: } else if (gpval<0 && xval>=xuptr[i]){
386: gpval = 0;
387: }
388: gpptr[i] = gpval;
390: /*
391: if (xlptr[i] >= xuptr[i]){
392: gpptr[i]=0.0;
393: } else if (xptr[i] <= xlptr[i]+eps){
394: gpptr[i] = PetscMin(gptr[i],zero);
395: } else if (xptr[i] >= xuptr[i]-eps){
396: gpptr[i] = PetscMax(gptr[i],zero);
397: } else {
398: gpptr[i] = gptr[i];
399: }
400: */
401: }
403: info=VecRestoreArray(X,&xptr); CHKERRQ(info);
404: info=VecRestoreArray(XL,&xlptr); CHKERRQ(info);
405: info=VecRestoreArray(XU,&xuptr); CHKERRQ(info);
406: info=VecRestoreArray(G,&gptr); CHKERRQ(info);
407: if (G!=GP){
408: info=VecRestoreArray(GP,&gpptr); CHKERRQ(info);
409: }
410: return(0);
411: }
415: int VecISAXPY(Vec Y, PetscScalar alpha, Vec X, IS YIS){
417: int i,info,in1,in2,xlow,xhigh,ylow,yhigh,*yis;
418: PetscScalar *x,*y;
421: info=ISGetLocalSize(YIS,&in1); CHKERRQ(info);
422: info=VecGetLocalSize(X,&in2); CHKERRQ(info);
424: if ( in1 != in2 )
425: SETERRQ(1,"Index set and X vector must be identically loaded over processors");
427: info=VecGetOwnershipRange(X, &xlow, &xhigh); CHKERRQ(info);
428: info=VecGetOwnershipRange(Y, &ylow, &yhigh); CHKERRQ(info);
430: info=ISGetIndices(YIS, &yis); CHKERRQ(info);
432: info=VecGetArray(X,&x); CHKERRQ(info);
433: if (X != Y){
434: info=VecGetArray(Y,&y); CHKERRQ(info);
435: } else {
436: y=x;
437: }
439: for (i=0; i<in1; i++){
440: if (yis[i]>=ylow && yis[i]<yhigh && i+xlow < xhigh){
441: y[yis[i]-ylow] = y[yis[i]-ylow] + (alpha)*x[i];
442: } else {
443: SETERRQ(1,"IS index out of range");
444: }
445: }
446:
447: info=ISRestoreIndices(YIS, &yis); CHKERRQ(info);
449: info=VecRestoreArray(X,&x); CHKERRQ(info);
450: if (X != Y){
451: info=VecRestoreArray(Y,&y); CHKERRQ(info);
452: }
454: info = PetscLogFlops(2*in1); CHKERRQ(info);
456: return(0);
457: }
461: int VecCreateSubVec(Vec A,IS S,Vec *B)
462: {
463: int info,nloc,n;
464: int i,low,high,nnloc;
465: int *ptrS;
466: PetscScalar zero=0.0;
467: PetscScalar *ptrB,*ptrA;
468: VecType type_name;
469: MPI_Comm comm;
474: info=ISGetLocalSize(S,&nloc);CHKERRQ(info);
475: info=ISGetSize(S,&n);CHKERRQ(info);
476: info=PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(info);
477: info=VecCreate(comm,B);CHKERRQ(info);
478: info=VecSetSizes(*B,nloc,n);CHKERRQ(info);
479: info=VecGetType(A,&type_name);CHKERRQ(info);
480: info=VecSetType(*B,type_name);CHKERRQ(info);
481: info=VecSet(*B, zero);CHKERRQ(info);
482:
483: info=VecGetOwnershipRange(A,&low,&high);CHKERRQ(info);
484: info=VecGetLocalSize(A,&nnloc);CHKERRQ(info);
485: info=VecGetArray(A,&ptrA);CHKERRQ(info);
486: info=VecGetArray(*B,&ptrB);CHKERRQ(info);
487: info=ISGetIndices(S,&ptrS);CHKERRQ(info);
488: for (i=0;i<nloc;i++){ ptrB[i]=ptrA[ptrS[i]-low]; }
489: info=VecRestoreArray(A,&ptrA);CHKERRQ(info);
490: info=VecRestoreArray(*B,&ptrB);CHKERRQ(info);
491: info=ISRestoreIndices(S,&ptrS);CHKERRQ(info);
493: return(0);
494: }
498: int VecStepMax(Vec X, Vec DX, PetscReal*step){
499: int i,info,nn;
500: PetscScalar stepmax=1.0e300;
501: PetscScalar *xx,*dx;
502: double t1,t2;
505: info = VecGetLocalSize(X,&nn);CHKERRQ(info);
506: info = VecGetArray(X,&xx);CHKERRQ(info);
507: info = VecGetArray(DX,&dx);CHKERRQ(info);
508: for (i=0;i<nn;i++){
509: if (xx[i] < 0){
510: SETERRQ(1,"Vector must be positive");
511: } else if (dx[i]<0){ stepmax=PetscMin(stepmax,-xx[i]/dx[i]);
512: }
513: }
514: info = VecRestoreArray(X,&xx);CHKERRQ(info);
515: info = VecRestoreArray(DX,&dx);CHKERRQ(info);
516: t1=(double)stepmax;
517: info = MPI_Allreduce(&t1,&t2,1,MPI_DOUBLE,MPI_MIN,X->comm);
518: CHKERRQ(info);
519: *step=(PetscReal)t2;
520: return(0);
521: }
525: int MatCreateProjectionOperator(Vec A, IS S , Vec B, Mat *MM){
527: int m,n,M,N,low,high;
528: int info,i,*iptr;
529: PetscScalar one=1.0;
530: Mat MMM;
538: info = VecGetSize(B,&M);CHKERRQ(info);
539: info = VecGetLocalSize(B,&m);CHKERRQ(info);
540: info = VecGetSize(A,&N);CHKERRQ(info);
541: info = VecGetLocalSize(A,&n);CHKERRQ(info);
543: info = MatCreateMPIAIJ(A->comm,m,n,M,N,1,PETSC_NULL,1,PETSC_NULL,&MMM);
544: CHKERRQ(info);
545: /*
546: info = MatCreateMPIBDiag(pv->comm,m,M,N,M,1,PETSC_NULL,PETSC_NULL,&A);
547: */
548:
549: info = ISGetSize(S,&n);CHKERRQ(info);
550: if (n!=M){
551: SETERRQ(1,"Size of index set does not match the second vector.");
552: }
553:
554: info = ISGetLocalSize(S,&n);CHKERRQ(info);
555: if (n!=m){
556: SETERRQ(1,"Local size of index set does not match the second vector.");
557: }
558: info=VecGetOwnershipRange(B,&low,&high);CHKERRQ(info);
559: info = ISGetIndices(S,&iptr);CHKERRQ(info);
560: for (i=0; i<n; i++){
561: info=MatSetValue(MMM,low+i,iptr[i],one,INSERT_VALUES);CHKERRQ(info);
562: }
563: info = ISRestoreIndices(S,&iptr);CHKERRQ(info);
565: info = MatAssemblyBegin(MMM,MAT_FINAL_ASSEMBLY);CHKERRQ(info);
566: info = MatAssemblyEnd(MMM,MAT_FINAL_ASSEMBLY);CHKERRQ(info);
568: info = MatCompress(MMM);CHKERRQ(info);
570: *MM = MMM;
572: return(0);
573: }
577: int VecStepMax2(Vec X, Vec DX, Vec XL, Vec XU, PetscReal *step2){
579: int i,nn1,nn2,info;
580: PetscScalar *xx,*dx,*xl,*xu;
581: PetscScalar stepmax2=0;
582: double t1,t2;
584: TaoFunctionBegin;
585: info = VecGetArray(X,&xx);CHKERRQ(info);
586: info = VecGetArray(XL,&xl);CHKERRQ(info);
587: info = VecGetArray(XU,&xu);CHKERRQ(info);
588: info = VecGetArray(DX,&dx);CHKERRQ(info);
589: info = VecGetLocalSize(X,&nn1);CHKERRQ(info);
590: info = VecGetLocalSize(DX,&nn2);CHKERRQ(info);
592: for (i=0;i<nn1;i++){
593: if (dx[i] > 0){
594: stepmax2=PetscMax(stepmax2,(xu[i]-xx[i])/dx[i]);
595: } else if (dx[i]<0){
596: stepmax2=PetscMax(stepmax2,(xl[i]-xx[i])/dx[i]);
597: }
598: }
599: info = VecRestoreArray(X,&xx);CHKERRQ(info);
600: info = VecRestoreArray(XL,&xl);CHKERRQ(info);
601: info = VecRestoreArray(XU,&xu);CHKERRQ(info);
602: info = VecRestoreArray(DX,&dx);CHKERRQ(info);
604: t1=(double)stepmax2;
605: info = MPI_Allreduce(&t1,&t2,1,MPI_DOUBLE,MPI_MAX,X->comm);
606: CHKERRQ(info);
607: *step2=(PetscReal)t2;
609: TaoFunctionReturn(0);
610: }
615: int VecBoundStepInfo(Vec X, Vec XL, Vec XU, Vec DX, PetscReal *boundmin, PetscReal *wolfemin, PetscReal *boundmax){
616: int info,i,n;
617: PetscScalar *x,*xl,*xu,*dx;
618: double tt,t1=1.0e+20, t2=1.0e+20, t3=0;
619: double tt1,tt2,tt3;
620: MPI_Comm comm;
621:
623: info=VecGetArray(X,&x);CHKERRQ(info);
624: info=VecGetArray(XL,&xl);CHKERRQ(info);
625: info=VecGetArray(XU,&xu);CHKERRQ(info);
626: info=VecGetArray(DX,&dx);CHKERRQ(info);
627: info = VecGetLocalSize(X,&n);CHKERRQ(info);
628: for (i=0;i<n;i++){
629: if (dx[i]>0){
630: tt=(xu[i]-x[i])/dx[i];
631: t1=PetscMin(t1,tt);
632: if (t1>0){
633: t2=PetscMin(t2,tt);
634: }
635: t3=PetscMax(t3,tt);
636: } else if (dx[i]<0){
637: tt=(xl[i]-x[i])/dx[i];
638: t1=PetscMin(t1,tt);
639: if (t1>0){
640: t2=PetscMin(t2,tt);
641: }
642: t3=PetscMax(t3,tt);
643: }
644: }
645: info=VecRestoreArray(X,&x);CHKERRQ(info);
646: info=VecRestoreArray(XL,&xl);CHKERRQ(info);
647: info=VecRestoreArray(XU,&xu);CHKERRQ(info);
648: info=VecRestoreArray(DX,&dx);CHKERRQ(info);
649: info=PetscObjectGetComm((PetscObject)X,&comm);CHKERRQ(info);
650:
651: if (boundmin){ info = MPI_Allreduce(&t1,&tt1,1,MPI_DOUBLE,MPI_MIN,comm);CHKERRQ(info);}
652: if (wolfemin){ info = MPI_Allreduce(&t2,&tt2,1,MPI_DOUBLE,MPI_MIN,comm);CHKERRQ(info);}
653: if (boundmax) { info = MPI_Allreduce(&t3,&tt3,1,MPI_DOUBLE,MPI_MAX,comm);CHKERRQ(info);}
655: *boundmin=(PetscReal)tt1;
656: *wolfemin=(PetscReal)tt2;
657: *boundmax=(PetscReal)tt3;
658: info = PetscInfo3(X,"Step Bound Info: Closest Bound: %6.4e, Wolfe: %6.4e, Max: %6.4e \n",*boundmin,*wolfemin,*boundmax); CHKERRQ(info);
659: return(0);
660: }
664: /*@C
665: VecISSetToConstant - Sets the elements of a vector, specified by an index set, to a constant
667: Input Parameter:
668: + S - a PETSc IS
669: . c - the constant
670: - V - a Vec
672: .seealso VecSet()
674: Level: advanced
675: @*/
676: int VecISSetToConstant(IS S, PetscScalar c, Vec V){
677: int info,i,nloc,low,high,*s;
678: PetscScalar *v;
685: info = VecGetOwnershipRange(V, &low, &high); CHKERRQ(info);
686: info = ISGetLocalSize(S,&nloc);CHKERRQ(info);
688: info = ISGetIndices(S, &s); CHKERRQ(info);
689: info = VecGetArray(V,&v); CHKERRQ(info);
690: for (i=0; i<nloc; i++){
691: v[s[i]-low] = c;
692: /*
693: if (s[i]>=low && s[i]<high){
694: v[s[i]-low] = c;
695: } else {
696: SETERRQ(1,"IS index out of range");
697: }
698: */
699: }
700:
701: info = ISRestoreIndices(S, &s); CHKERRQ(info);
702: info = VecRestoreArray(V,&v); CHKERRQ(info);
704: return(0);
706: }
710: int VecPow(Vec Vec1, double p)
711: {
712: int i,n,info;
713: PetscScalar *v1;
718: info = VecGetArray(Vec1, &v1); CHKERRQ(info);
719: info = VecGetLocalSize(Vec1, &n); CHKERRQ(info);
721: if (1.0 == p) {
722: }
723: else if (-1.0 == p) {
724: for (i = 0; i < n; ++i){
725: v1[i] = 1.0 / v1[i];
726: }
727: }
728: else if (0.0 == p) {
729: for (i = 0; i < n; ++i){
730: // Not-a-number left alone
731: // Infinity set to one
732: if (v1[i] == v1[i]) {
733: v1[i] = 1.0;
734: }
735: }
736: }
737: else if (0.5 == p) {
738: for (i = 0; i < n; ++i) {
739: if (v1[i] >= 0) {
740: v1[i] = sqrt(v1[i]);
741: }
742: else {
743: v1[i] = TAO_INFINITY;
744: }
745: }
746: }
747: else if (-0.5 == p) {
748: for (i = 0; i < n; ++i) {
749: if (v1[i] >= 0) {
750: v1[i] = 1.0 / sqrt(v1[i]);
751: }
752: else {
753: v1[i] = TAO_INFINITY;
754: }
755: }
756: }
757: else if (2.0 == p) {
758: for (i = 0; i < n; ++i){
759: v1[i] *= v1[i];
760: }
761: }
762: else if (-2.0 == p) {
763: for (i = 0; i < n; ++i){
764: v1[i] = 1.0 / (v1[i] * v1[i]);
765: }
766: }
767: else {
768: for (i = 0; i < n; ++i) {
769: if (v1[i] >= 0) {
770: v1[i] = pow(v1[i], p);
771: }
772: else {
773: v1[i] = TAO_INFINITY;
774: }
775: }
776: }
778: info = VecRestoreArray(Vec1,&v1); CHKERRQ(info);
779: return(0);
780: }