Actual source code: bvec2.c

petsc-3.4.2 2013-07-02
  2: /*
  3:    Implements the sequential vectors.
  4: */

  6: #include <../src/vec/vec/impls/dvecimpl.h>          /*I "petscvec.h" I*/
  7: #include <../src/vec/vec/impls/mpi/pvecimpl.h> /* For VecView_MPI_HDF5 */
  8: #include <petscblaslapack.h>
  9: #include <petscthreadcomm.h>

 11: #if defined(PETSC_HAVE_HDF5)
 12: extern PetscErrorCode VecView_MPI_HDF5(Vec,PetscViewer);
 13: #endif

 15: #if defined(PETSC_THREADCOMM_ACTIVE)
 16: PetscErrorCode VecPointwiseMax_kernel(PetscInt thread_id,Vec win,Vec xin,Vec yin)
 17: {
 18:   PetscErrorCode    ierr;
 19:   PetscInt          *trstarts=win->map->trstarts;
 20:   PetscInt          i;
 21:   const PetscScalar *xx,*yy;
 22:   PetscScalar       *ww;

 24:   VecGetArrayRead(xin,&xx);
 25:   VecGetArrayRead(yin,&yy);
 26:   VecGetArray(win,&ww);

 28:   for (i=trstarts[thread_id]; i < trstarts[thread_id+1];i++) ww[i] = PetscMax(PetscRealPart(xx[i]),PetscRealPart(yy[i]));

 30:   VecRestoreArrayRead(xin,&xx);
 31:   VecRestoreArrayRead(yin,&yy);
 32:   VecRestoreArray(win,&ww);
 33:   return 0;
 34: }

 38: PetscErrorCode VecPointwiseMax_Seq(Vec win,Vec xin,Vec yin)
 39: {

 43:   PetscThreadCommRunKernel3(PetscObjectComm((PetscObject)win),(PetscThreadKernel)VecPointwiseMax_kernel,win,xin,yin);
 44:   PetscLogFlops(win->map->n);
 45:   return(0);
 46: }
 47: #else
 50: PetscErrorCode VecPointwiseMax_Seq(Vec win,Vec xin,Vec yin)
 51: {
 53:   PetscInt       n = win->map->n,i;
 54:   PetscScalar    *ww,*xx,*yy; /* cannot make xx or yy const since might be ww */

 57:   VecGetArrayRead(xin,(const PetscScalar**)&xx);
 58:   VecGetArrayRead(yin,(const PetscScalar**)&yy);
 59:   VecGetArray(win,&ww);

 61:   for (i=0; i<n; i++) ww[i] = PetscMax(PetscRealPart(xx[i]),PetscRealPart(yy[i]));

 63:   VecRestoreArrayRead(xin,(const PetscScalar**)&xx);
 64:   VecRestoreArrayRead(yin,(const PetscScalar**)&yy);
 65:   VecRestoreArray(win,&ww);
 66:   PetscLogFlops(n);
 67:   return(0);
 68: }
 69: #endif

 71: #if defined(PETSC_THREADCOMM_ACTIVE)
 72: PetscErrorCode VecPointwiseMin_kernel(PetscInt thread_id,Vec win,Vec xin,Vec yin)
 73: {
 74:   PetscErrorCode    ierr;
 75:   PetscInt          *trstarts=win->map->trstarts;
 76:   PetscInt          i;
 77:   const PetscScalar *xx,*yy;
 78:   PetscScalar       *ww;

 80:   VecGetArrayRead(xin,&xx);
 81:   VecGetArrayRead(yin,&yy);
 82:   VecGetArray(win,&ww);

 84:   for (i=trstarts[thread_id]; i < trstarts[thread_id+1];i++) ww[i] = PetscMin(PetscRealPart(xx[i]),PetscRealPart(yy[i]));

 86:   VecRestoreArrayRead(xin,&xx);
 87:   VecRestoreArrayRead(yin,&yy);
 88:   VecRestoreArray(win,&ww);
 89:   return 0;
 90: }

 94: PetscErrorCode VecPointwiseMin_Seq(Vec win,Vec xin,Vec yin)
 95: {

 99:   PetscThreadCommRunKernel3(PetscObjectComm((PetscObject)win),(PetscThreadKernel)VecPointwiseMin_kernel,win,xin,yin);
100:   PetscLogFlops(win->map->n);
101:   return(0);
102: }
103: #else
106: PetscErrorCode VecPointwiseMin_Seq(Vec win,Vec xin,Vec yin)
107: {
109:   PetscInt       n = win->map->n,i;
110:   PetscScalar    *ww,*xx,*yy; /* cannot make xx or yy const since might be ww */

113:   VecGetArrayRead(xin,(const PetscScalar**)&xx);
114:   VecGetArrayRead(yin,(const PetscScalar**)&yy);
115:   VecGetArray(win,&ww);

117:   for (i=0; i<n; i++) ww[i] = PetscMin(PetscRealPart(xx[i]),PetscRealPart(yy[i]));

119:   VecRestoreArrayRead(xin,(const PetscScalar**)&xx);
120:   VecRestoreArrayRead(yin,(const PetscScalar**)&yy);
121:   VecRestoreArray(win,&ww);
122:   PetscLogFlops(n);
123:   return(0);
124: }
125: #endif

127: #if defined(PETSC_THREADCOMM_ACTIVE)
128: PetscErrorCode VecPointwiseMaxAbs_kernel(PetscInt thread_id,Vec win,Vec xin,Vec yin)
129: {
130:   PetscErrorCode    ierr;
131:   PetscInt          *trstarts=win->map->trstarts;
132:   PetscInt          i;
133:   const PetscScalar *xx,*yy;
134:   PetscScalar       *ww;

136:   VecGetArrayRead(xin,&xx);
137:   VecGetArrayRead(yin,&yy);
138:   VecGetArray(win,&ww);

140:   for (i=trstarts[thread_id]; i < trstarts[thread_id+1];i++) ww[i] = PetscMax(PetscAbsScalar(xx[i]),PetscAbsScalar(yy[i]));

142:   VecRestoreArrayRead(xin,&xx);
143:   VecRestoreArrayRead(yin,&yy);
144:   VecRestoreArray(win,&ww);
145:   return 0;
146: }

150: PetscErrorCode VecPointwiseMaxAbs_Seq(Vec win,Vec xin,Vec yin)
151: {

155:   PetscThreadCommRunKernel3(PetscObjectComm((PetscObject)win),(PetscThreadKernel)VecPointwiseMaxAbs_kernel,win,xin,yin);
156:   PetscLogFlops(win->map->n);
157:   return(0);
158: }
159: #else
162: PetscErrorCode VecPointwiseMaxAbs_Seq(Vec win,Vec xin,Vec yin)
163: {
165:   PetscInt       n = win->map->n,i;
166:   PetscScalar    *ww,*xx,*yy; /* cannot make xx or yy const since might be ww */

169:   VecGetArrayRead(xin,(const PetscScalar**)&xx);
170:   VecGetArrayRead(yin,(const PetscScalar**)&yy);
171:   VecGetArray(win,&ww);

173:   for (i=0; i<n; i++) ww[i] = PetscMax(PetscAbsScalar(xx[i]),PetscAbsScalar(yy[i]));

175:   PetscLogFlops(n);
176:   VecRestoreArrayRead(xin,(const PetscScalar**)&xx);
177:   VecRestoreArrayRead(yin,(const PetscScalar**)&yy);
178:   VecRestoreArray(win,&ww);
179:   return(0);
180: }
181: #endif

183: #include <../src/vec/vec/impls/seq/ftn-kernels/fxtimesy.h>
184: #if defined(PETSC_THREADCOMM_ACTIVE)
185: PetscErrorCode VecPointwiseMult_kernel(PetscInt thread_id,Vec win,Vec xin,Vec yin)
186: {
187:   PetscErrorCode    ierr;
188:   PetscInt          *trstarts=win->map->trstarts;
189:   PetscInt          i;
190:   const PetscScalar *xx,*yy;
191:   PetscScalar       *ww;

193:   VecGetArrayRead(xin,&xx);
194:   VecGetArrayRead(yin,&yy);
195:   VecGetArray(win,&ww);
196:   if (ww == xx) {
197:     for (i=trstarts[thread_id]; i < trstarts[thread_id+1]; i++) ww[i] *= yy[i];
198:   } else if (ww == yy) {
199:     for (i=trstarts[thread_id]; i < trstarts[thread_id+1]; i++) ww[i] *= xx[i];
200:   } else {
201: #if defined(PETSC_USE_FORTRAN_KERNEL_XTIMESY)
202:     PetscInt start,n;
203:     start = trstarts[thread_id];
204:     n     = trstarts[thread_id+1] - trstarts[thread_id];
205:     fortranxtimesy_(xx+start,yy+start,ww+start,&n);
206:   }
207: #else
208:     for (i=trstarts[thread_id]; i < trstarts[thread_id+1]; i++) ww[i] = xx[i] * yy[i];
209:   }
210: #endif
211:   VecRestoreArrayRead(xin,&xx);
212:   VecRestoreArrayRead(yin,&yy);
213:   VecRestoreArray(win,&ww);
214:   return 0;
215: }

219: PetscErrorCode VecPointwiseMult_Seq(Vec win,Vec xin,Vec yin)
220: {

224:   PetscThreadCommRunKernel3(PetscObjectComm((PetscObject)win),(PetscThreadKernel)VecPointwiseMult_kernel,win,xin,yin);
225:   PetscLogFlops(win->map->n);
226:   return(0);
227: }
228: #else
231: PetscErrorCode VecPointwiseMult_Seq(Vec win,Vec xin,Vec yin)
232: {
234:   PetscInt       n = win->map->n,i;
235:   PetscScalar    *ww,*xx,*yy; /* cannot make xx or yy const since might be ww */

238:   VecGetArrayRead(xin,(const PetscScalar**)&xx);
239:   VecGetArrayRead(yin,(const PetscScalar**)&yy);
240:   VecGetArray(win,&ww);
241:   if (ww == xx) {
242:     for (i=0; i<n; i++) ww[i] *= yy[i];
243:   } else if (ww == yy) {
244:     for (i=0; i<n; i++) ww[i] *= xx[i];
245:   } else {
246: #if defined(PETSC_USE_FORTRAN_KERNEL_XTIMESY)
247:     fortranxtimesy_(xx,yy,ww,&n);
248: #else
249:     for (i=0; i<n; i++) ww[i] = xx[i] * yy[i];
250: #endif
251:   }
252:   VecRestoreArrayRead(xin,(const PetscScalar**)&xx);
253:   VecRestoreArrayRead(yin,(const PetscScalar**)&yy);
254:   VecRestoreArray(win,&ww);
255:   PetscLogFlops(n);
256:   return(0);
257: }
258: #endif

260: #if defined(PETSC_THREADCOMM_ACTIVE)
261: PetscErrorCode VecPointwiseDivide_kernel(PetscInt thread_id,Vec win,Vec xin,Vec yin)
262: {
263:   PetscErrorCode    ierr;
264:   PetscInt          *trstarts=win->map->trstarts;
265:   PetscInt          i;
266:   const PetscScalar *xx,*yy;
267:   PetscScalar       *ww;

269:   VecGetArrayRead(xin,&xx);
270:   VecGetArrayRead(yin,&yy);
271:   VecGetArray(win,&ww);

273:   for (i=trstarts[thread_id]; i < trstarts[thread_id+1]; i++) ww[i] = xx[i] / yy[i];

275:   VecRestoreArrayRead(xin,&xx);
276:   VecRestoreArrayRead(yin,&yy);
277:   VecRestoreArray(win,&ww);
278:   return 0;
279: }

283: PetscErrorCode VecPointwiseDivide_Seq(Vec win,Vec xin,Vec yin)
284: {

288:   PetscThreadCommRunKernel3(PetscObjectComm((PetscObject)win),(PetscThreadKernel)VecPointwiseDivide_kernel,win,xin,yin);
289:   PetscLogFlops(win->map->n);
290:   return(0);
291: }
292: #else
295: PetscErrorCode VecPointwiseDivide_Seq(Vec win,Vec xin,Vec yin)
296: {
298:   PetscInt       n = win->map->n,i;
299:   PetscScalar    *ww,*xx,*yy; /* cannot make xx or yy const since might be ww */

302:   VecGetArrayRead(xin,(const PetscScalar**)&xx);
303:   VecGetArrayRead(yin,(const PetscScalar**)&yy);
304:   VecGetArray(win,&ww);

306:   for (i=0; i<n; i++) ww[i] = xx[i] / yy[i];

308:   PetscLogFlops(n);
309:   VecRestoreArrayRead(xin,(const PetscScalar**)&xx);
310:   VecRestoreArrayRead(yin,(const PetscScalar**)&yy);
311:   VecRestoreArray(win,&ww);
312:   return(0);
313: }
314: #endif

316: #if defined(PETSC_THREADCOMM_ACTIVE)
317: PetscErrorCode VecSetRandom_kernel(PetscInt thread_id,Vec xin, PetscRandom r)
318: {
320:   PetscInt       *trstarts=xin->map->trstarts;
321:   PetscInt       i;
322:   PetscScalar    *xx;

324:   VecGetArray(xin,&xx);
325:   for (i=trstarts[thread_id]; i < trstarts[thread_id+1]; i++) {
326:     PetscRandomGetValue(r,&xx[i]);
327:   }
328:   VecRestoreArray(xin,&xx);
329:   return 0;
330: }

334: PetscErrorCode VecSetRandom_Seq(Vec xin,PetscRandom r)
335: {

339:   PetscThreadCommRunKernel2(PetscObjectComm((PetscObject)xin),(PetscThreadKernel)VecSetRandom_kernel,xin,r);
340:   return(0);
341: }
342: #else
345: PetscErrorCode VecSetRandom_Seq(Vec xin,PetscRandom r)
346: {
348:   PetscInt       n = xin->map->n,i;
349:   PetscScalar    *xx;

352:   VecGetArray(xin,&xx);
353:   for (i=0; i<n; i++) {PetscRandomGetValue(r,&xx[i]);}
354:   VecRestoreArray(xin,&xx);
355:   return(0);
356: }
357: #endif

361: PetscErrorCode VecGetSize_Seq(Vec vin,PetscInt *size)
362: {
364:   *size = vin->map->n;
365:   return(0);
366: }


369: #if defined(PETSC_THREADCOMM_ACTIVE)
370: PetscErrorCode VecConjugate_kernel(PetscInt thread_id,Vec xin)
371: {
373:   PetscScalar    *x;
374:   PetscInt       *trstarts=xin->map->trstarts,i;

376:   VecGetArray(xin,&x);
377:   for (i=trstarts[thread_id]; i < trstarts[thread_id+1]; i++) x[i] = PetscConj(x[i]);
378:   VecRestoreArray(xin,&x);
379:   return 0;
380: }

384: PetscErrorCode VecConjugate_Seq(Vec xin)
385: {

389:   PetscThreadCommRunKernel1(PetscObjectComm((PetscObject)xin),(PetscThreadKernel)VecConjugate_kernel,xin);
390:   return(0);
391: }
392: #else
395: PetscErrorCode VecConjugate_Seq(Vec xin)
396: {
397:   PetscScalar    *x;
398:   PetscInt       n = xin->map->n;

402:   VecGetArray(xin,&x);
403:   while (n-->0) {
404:     *x = PetscConj(*x);
405:     x++;
406:   }
407:   VecRestoreArray(xin,&x);
408:   return(0);
409: }
410: #endif

414: PetscErrorCode VecResetArray_Seq(Vec vin)
415: {
416:   Vec_Seq *v = (Vec_Seq*)vin->data;

419:   v->array         = v->unplacedarray;
420:   v->unplacedarray = 0;
421:   return(0);
422: }

424: #if defined(PETSC_THREADCOMM_ACTIVE)
425: PetscErrorCode VecCopy_kernel(PetscInt thread_id,Vec xin,Vec yin)
426: {
427:   PetscErrorCode    ierr;
428:   PetscInt          *trstarts=yin->map->trstarts;
429:   PetscScalar       *ya;
430:   const PetscScalar *xa;
431:   PetscInt          start,end,n;

433:   VecGetArrayRead(xin,&xa);
434:   VecGetArray(yin,&ya);
435:   start = trstarts[thread_id];
436:   end   = trstarts[thread_id+1];
437:   n     = end-start;
438:   PetscMemcpy(ya+start,xa+start,n*sizeof(PetscScalar));
439:   VecRestoreArrayRead(xin,&xa);
440:   VecRestoreArray(yin,&ya);
441:   return 0;
442: }

446: PetscErrorCode VecCopy_Seq(Vec xin,Vec yin)
447: {

451:   if (xin != yin) {
452:     PetscThreadCommRunKernel2(PetscObjectComm((PetscObject)xin),(PetscThreadKernel)VecCopy_kernel,xin,yin);
453:   }
454:   return(0);
455: }
456: #else
459: PetscErrorCode VecCopy_Seq(Vec xin,Vec yin)
460: {
461:   PetscScalar       *ya;
462:   const PetscScalar *xa;
463:   PetscErrorCode    ierr;

466:   if (xin != yin) {
467:     VecGetArrayRead(xin,&xa);
468:     VecGetArray(yin,&ya);
469:     PetscMemcpy(ya,xa,xin->map->n*sizeof(PetscScalar));
470:     VecRestoreArrayRead(xin,&xa);
471:     VecRestoreArray(yin,&ya);
472:   }
473:   return(0);
474: }
475: #endif

477: #if defined(PETSC_THREADCOMM_ACTIVE)
478: PetscErrorCode VecSwap_kernel(PetscInt thread_id,Vec xin,Vec yin)
479: {
480:   PetscErrorCode   ierr;
481:   PetscInt         *trstarts=xin->map->trstarts;
482:   PetscInt         start,end,n;
483:   PetscBLASInt     one=1,bn;
484:   PetscScalar      *xa,*ya;

486:   VecGetArray(xin,&xa);
487:   VecGetArray(yin,&ya);
488:   start = trstarts[thread_id];
489:   end   = trstarts[thread_id+1];
490:   n     = end-start;
491:   PetscBLASIntCast(n,&bn);
492:   PetscStackCallBLAS("BLASswap",BLASswap_(&bn,xa+start,&one,ya+start,&one));
493:   VecRestoreArray(xin,&xa);
494:   VecRestoreArray(yin,&ya);
495:   return 0;
496: }

500: PetscErrorCode VecSwap_Seq(Vec xin,Vec yin)
501: {

505:   if (xin != yin) {
506:     PetscThreadCommRunKernel2(PetscObjectComm((PetscObject)xin),(PetscThreadKernel)VecSwap_kernel,xin,yin);
507:   }
508:   return(0);
509: }
510: #else

514: PetscErrorCode VecSwap_Seq(Vec xin,Vec yin)
515: {
516:   PetscScalar    *ya, *xa;
518:   PetscBLASInt   one = 1,bn;

521:   if (xin != yin) {
522:     PetscBLASIntCast(xin->map->n,&bn);
523:     VecGetArray(xin,&xa);
524:     VecGetArray(yin,&ya);
525:     PetscStackCallBLAS("BLASswap",BLASswap_(&bn,xa,&one,ya,&one));
526:     VecRestoreArray(xin,&xa);
527:     VecRestoreArray(yin,&ya);
528:   }
529:   return(0);
530: }
531: #endif

533: #include <../src/vec/vec/impls/seq/ftn-kernels/fnorm.h>
534: #if defined(PETSC_THREADCOMM_ACTIVE)
535: PetscErrorCode VecNorm_kernel(PetscInt thread_id,Vec xin,NormType* type_p,PetscThreadCommReduction red)
536: {
538:   PetscInt       *trstarts=xin->map->trstarts;
539:   PetscInt       start,end,n;
540:   PetscBLASInt   one = 1,bn;
541:   const PetscScalar *xx;
542:   PetscReal      z_loc;
543:   NormType       type=*type_p;

545:   start = trstarts[thread_id];
546:   end   = trstarts[thread_id+1];
547:   n     = end - start;
548:   PetscBLASIntCast(n,&bn);
549:   VecGetArrayRead(xin,&xx);
550:   if (type == NORM_2 || type == NORM_FROBENIUS) {
551:     z_loc = PetscRealPart(BLASdot_(&bn,xx+start,&one,xx+start,&one));
552:     PetscThreadReductionKernelPost(thread_id,red,(void*)&z_loc);
553:   } else if (type == NORM_INFINITY) {
554:     PetscInt i;
555:     PetscReal max=0.0,tmp;
556:     for (i=trstarts[thread_id]; i < trstarts[thread_id+1]; i++) {
557:       if ((tmp = PetscAbsScalar(xx[i])) > max) max = tmp;
558:       /* check special case of tmp == NaN */
559:       if (tmp != tmp) {max = tmp; break;}
560:     }
561:     PetscThreadReductionKernelPost(thread_id,red,(void*)&max);
562:   } else if (type == NORM_1) {
563:     PetscStackCallBLAS("BLASasum",z_loc = BLASasum_(&bn,xx+start,&one));
564:     PetscThreadReductionKernelPost(thread_id,red,(void*)&z_loc);
565:   }
566:   VecRestoreArrayRead(xin,&xx);
567:   return 0;
568: }

572: PetscErrorCode VecNorm_Seq(Vec xin,NormType type,PetscReal *z)
573: {
574:   PetscErrorCode           ierr;
575:   PetscThreadCommReduction red;

578:   if (type == NORM_2 || type == NORM_FROBENIUS) {
579:     PetscThreadReductionBegin(PetscObjectComm((PetscObject)xin),THREADCOMM_SUM,PETSC_REAL,1,&red);
580:     PetscThreadCommRunKernel3(PetscObjectComm((PetscObject)xin),(PetscThreadKernel)VecNorm_kernel,xin,(void*)&type,red);
581:     PetscThreadReductionEnd(red,z);
582:     *z   = PetscSqrtReal(*z);
583:   } else if (type == NORM_INFINITY) {
584:     PetscThreadReductionBegin(PetscObjectComm((PetscObject)xin),THREADCOMM_MAX,PETSC_REAL,1,&red);
585:     PetscThreadCommRunKernel3(PetscObjectComm((PetscObject)xin),(PetscThreadKernel)VecNorm_kernel,xin,(void*)&type,red);
586:     PetscThreadReductionEnd(red,z);
587:   } else if (type == NORM_1) {
588:     PetscThreadReductionBegin(PetscObjectComm((PetscObject)xin),THREADCOMM_SUM,PETSC_REAL,1,&red);
589:     PetscThreadCommRunKernel3(PetscObjectComm((PetscObject)xin),(PetscThreadKernel)VecNorm_kernel,xin,(void*)&type,red);
590:     PetscThreadReductionEnd(red,z);
591:   } else if (type == NORM_1_AND_2) {
592:     VecNorm_Seq(xin,NORM_1,z);
593:     VecNorm_Seq(xin,NORM_2,z+1);
594:   }
595:   return(0);
596: }
597: #else

601: PetscErrorCode VecNorm_Seq(Vec xin,NormType type,PetscReal *z)
602: {
603:   const PetscScalar *xx;
604:   PetscErrorCode    ierr;
605:   PetscInt          n = xin->map->n;
606:   PetscBLASInt      one = 1, bn;

609:   PetscBLASIntCast(n,&bn);
610:   if (type == NORM_2 || type == NORM_FROBENIUS) {
611:     VecGetArrayRead(xin,&xx);
612:     *z   = PetscRealPart(BLASdot_(&bn,xx,&one,xx,&one));
613:     *z   = PetscSqrtReal(*z);
614:     VecRestoreArrayRead(xin,&xx);
615:     PetscLogFlops(PetscMax(2.0*n-1,0.0));
616:   } else if (type == NORM_INFINITY) {
617:     PetscInt  i;
618:     PetscReal max = 0.0,tmp;

620:     VecGetArrayRead(xin,&xx);
621:     for (i=0; i<n; i++) {
622:       if ((tmp = PetscAbsScalar(*xx)) > max) max = tmp;
623:       /* check special case of tmp == NaN */
624:       if (tmp != tmp) {max = tmp; break;}
625:       xx++;
626:     }
627:     VecRestoreArrayRead(xin,&xx);
628:     *z   = max;
629:   } else if (type == NORM_1) {
630:     VecGetArrayRead(xin,&xx);
631:     PetscStackCallBLAS("BLASasum",*z   = BLASasum_(&bn,xx,&one));
632:     VecRestoreArrayRead(xin,&xx);
633:     PetscLogFlops(PetscMax(n-1.0,0.0));
634:   } else if (type == NORM_1_AND_2) {
635:     VecNorm_Seq(xin,NORM_1,z);
636:     VecNorm_Seq(xin,NORM_2,z+1);
637:   }
638:   return(0);
639: }
640: #endif

644: PetscErrorCode VecView_Seq_ASCII(Vec xin,PetscViewer viewer)
645: {
646:   PetscErrorCode    ierr;
647:   PetscInt          i,n = xin->map->n;
648:   const char        *name;
649:   PetscViewerFormat format;
650:   const PetscScalar *xv;

653:   VecGetArrayRead(xin,&xv);
654:   PetscViewerGetFormat(viewer,&format);
655:   if (format == PETSC_VIEWER_ASCII_MATLAB) {
656:     PetscObjectGetName((PetscObject)xin,&name);
657:     PetscViewerASCIIPrintf(viewer,"%s = [\n",name);
658:     for (i=0; i<n; i++) {
659: #if defined(PETSC_USE_COMPLEX)
660:       if (PetscImaginaryPart(xv[i]) > 0.0) {
661:         PetscViewerASCIIPrintf(viewer,"%18.16e + %18.16ei\n",PetscRealPart(xv[i]),PetscImaginaryPart(xv[i]));
662:       } else if (PetscImaginaryPart(xv[i]) < 0.0) {
663:         PetscViewerASCIIPrintf(viewer,"%18.16e - %18.16ei\n",PetscRealPart(xv[i]),-PetscImaginaryPart(xv[i]));
664:       } else {
665:         PetscViewerASCIIPrintf(viewer,"%18.16e\n",PetscRealPart(xv[i]));
666:       }
667: #else
668:       PetscViewerASCIIPrintf(viewer,"%18.16e\n",(double) xv[i]);
669: #endif
670:     }
671:     PetscViewerASCIIPrintf(viewer,"];\n");
672:   } else if (format == PETSC_VIEWER_ASCII_SYMMODU) {
673:     for (i=0; i<n; i++) {
674: #if defined(PETSC_USE_COMPLEX)
675:       PetscViewerASCIIPrintf(viewer,"%18.16e %18.16e\n",PetscRealPart(xv[i]),PetscImaginaryPart(xv[i]));
676: #else
677:       PetscViewerASCIIPrintf(viewer,"%18.16e\n",(double)xv[i]);
678: #endif
679:     }
680:   } else if (format == PETSC_VIEWER_ASCII_VTK || format == PETSC_VIEWER_ASCII_VTK_CELL) {
681:     /*
682:        state 0: No header has been output
683:        state 1: Only POINT_DATA has been output
684:        state 2: Only CELL_DATA has been output
685:        state 3: Output both, POINT_DATA last
686:        state 4: Output both, CELL_DATA last
687:     */
688:     static PetscInt stateId = -1;
689:     int outputState = 0;
690:     PetscBool  hasState;
691:     int doOutput = 0;
692:     PetscInt bs, b;

694:     if (stateId < 0) {
695:       PetscObjectComposedDataRegister(&stateId);
696:     }
697:     PetscObjectComposedDataGetInt((PetscObject) viewer, stateId, outputState, hasState);
698:     if (!hasState) outputState = 0;
699:     PetscObjectGetName((PetscObject) xin, &name);
700:     VecGetBlockSize(xin, &bs);
701:     if ((bs < 1) || (bs > 3)) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE, "VTK can only handle 3D objects, but vector dimension is %d", bs);
702:     if (format == PETSC_VIEWER_ASCII_VTK) {
703:       if (outputState == 0) {
704:         outputState = 1;
705:         doOutput = 1;
706:       } else if (outputState == 1) doOutput = 0;
707:       else if (outputState == 2) {
708:         outputState = 3;
709:         doOutput = 1;
710:       } else if (outputState == 3) doOutput = 0;
711:       else if (outputState == 4) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE, "Tried to output POINT_DATA again after intervening CELL_DATA");

713:       if (doOutput) {
714:         PetscViewerASCIIPrintf(viewer, "POINT_DATA %d\n", n/bs);
715:       }
716:     } else {
717:       if (outputState == 0) {
718:         outputState = 2;
719:         doOutput = 1;
720:       } else if (outputState == 1) {
721:         outputState = 4;
722:         doOutput = 1;
723:       } else if (outputState == 2) doOutput = 0;
724:       else if (outputState == 3) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE, "Tried to output CELL_DATA again after intervening POINT_DATA");
725:       else if (outputState == 4) doOutput = 0;

727:       if (doOutput) {
728:         PetscViewerASCIIPrintf(viewer, "CELL_DATA %d\n", n);
729:       }
730:     }
731:     PetscObjectComposedDataSetInt((PetscObject) viewer, stateId, outputState);
732:     if (name) {
733:       if (bs == 3) {
734:         PetscViewerASCIIPrintf(viewer, "VECTORS %s double\n", name);
735:       } else {
736:         PetscViewerASCIIPrintf(viewer, "SCALARS %s double %d\n", name, bs);
737:       }
738:     } else {
739:       PetscViewerASCIIPrintf(viewer, "SCALARS scalars double %d\n", bs);
740:     }
741:     if (bs != 3) {
742:       PetscViewerASCIIPrintf(viewer, "LOOKUP_TABLE default\n");
743:     }
744:     for (i=0; i<n/bs; i++) {
745:       for (b=0; b<bs; b++) {
746:         if (b > 0) {
747:           PetscViewerASCIIPrintf(viewer," ");
748:         }
749: #if !defined(PETSC_USE_COMPLEX)
750:         PetscViewerASCIIPrintf(viewer,"%g",xv[i*bs+b]);
751: #endif
752:       }
753:       PetscViewerASCIIPrintf(viewer,"\n");
754:     }
755:   } else if (format == PETSC_VIEWER_ASCII_VTK_COORDS) {
756:     PetscInt bs, b;

758:     VecGetBlockSize(xin, &bs);
759:     if ((bs < 1) || (bs > 3)) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE, "VTK can only handle 3D objects, but vector dimension is %d", bs);
760:     for (i=0; i<n/bs; i++) {
761:       for (b=0; b<bs; b++) {
762:         if (b > 0) {
763:           PetscViewerASCIIPrintf(viewer," ");
764:         }
765: #if !defined(PETSC_USE_COMPLEX)
766:         PetscViewerASCIIPrintf(viewer,"%g",xv[i*bs+b]);
767: #endif
768:       }
769:       for (b=bs; b<3; b++) {
770:         PetscViewerASCIIPrintf(viewer," 0.0");
771:       }
772:       PetscViewerASCIIPrintf(viewer,"\n");
773:     }
774:   } else if (format == PETSC_VIEWER_ASCII_PCICE) {
775:     PetscInt bs, b;

777:     VecGetBlockSize(xin, &bs);
778:     if ((bs < 1) || (bs > 3)) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE, "PCICE can only handle up to 3D objects, but vector dimension is %d", bs);
779:     PetscViewerASCIIPrintf(viewer,"%D\n", xin->map->N/bs);
780:     for (i=0; i<n/bs; i++) {
781:       PetscViewerASCIIPrintf(viewer,"%7D   ", i+1);
782:       for (b=0; b<bs; b++) {
783:         if (b > 0) {
784:           PetscViewerASCIIPrintf(viewer," ");
785:         }
786: #if !defined(PETSC_USE_COMPLEX)
787:         PetscViewerASCIIPrintf(viewer,"% 12.5E",xv[i*bs+b]);
788: #endif
789:       }
790:       PetscViewerASCIIPrintf(viewer,"\n");
791:     }
792:   } else if (format == PETSC_VIEWER_ASCII_INFO || format == PETSC_VIEWER_ASCII_INFO_DETAIL) return(0);
793:   else {
794:     PetscObjectPrintClassNamePrefixType((PetscObject)xin,viewer,"Vector Object");
795:     for (i=0; i<n; i++) {
796:       if (format == PETSC_VIEWER_ASCII_INDEX) {
797:         PetscViewerASCIIPrintf(viewer,"%D: ",i);
798:       }
799: #if defined(PETSC_USE_COMPLEX)
800:       if (PetscImaginaryPart(xv[i]) > 0.0) {
801:         PetscViewerASCIIPrintf(viewer,"%g + %g i\n",PetscRealPart(xv[i]),PetscImaginaryPart(xv[i]));
802:       } else if (PetscImaginaryPart(xv[i]) < 0.0) {
803:         PetscViewerASCIIPrintf(viewer,"%g - %g i\n",PetscRealPart(xv[i]),-PetscImaginaryPart(xv[i]));
804:       } else {
805:         PetscViewerASCIIPrintf(viewer,"%g\n",PetscRealPart(xv[i]));
806:       }
807: #else
808:       PetscViewerASCIIPrintf(viewer,"%g\n",(double)xv[i]);
809: #endif
810:     }
811:   }
812:   PetscViewerFlush(viewer);
813:   VecRestoreArrayRead(xin,&xv);
814:   return(0);
815: }

817: #include <petscdraw.h>
820: PetscErrorCode VecView_Seq_Draw_LG(Vec xin,PetscViewer v)
821: {
822:   PetscErrorCode    ierr;
823:   PetscInt          i,c,bs = xin->map->bs,n = xin->map->n/bs;
824:   PetscDraw         win;
825:   PetscReal         *xx;
826:   PetscDrawLG       lg;
827:   const PetscScalar *xv;
828:   PetscReal         *yy;

831:   PetscMalloc(n*sizeof(PetscReal),&xx);
832:   PetscMalloc(n*sizeof(PetscReal),&yy);
833:   VecGetArrayRead(xin,&xv);
834:   for (c=0; c<bs; c++) {
835:     PetscViewerDrawGetDrawLG(v,c,&lg);
836:     PetscDrawLGGetDraw(lg,&win);
837:     PetscDrawCheckResizedWindow(win);
838:     PetscDrawLGReset(lg);
839:     for (i=0; i<n; i++) {
840:       xx[i] = (PetscReal) i;
841:       yy[i] = PetscRealPart(xv[c + i*bs]);
842:     }
843:     PetscDrawLGAddPoints(lg,n,&xx,&yy);
844:     PetscDrawLGDraw(lg);
845:     PetscDrawSynchronizedFlush(win);
846:   }
847:   VecRestoreArrayRead(xin,&xv);
848:   PetscFree(yy);
849:   PetscFree(xx);
850:   return(0);
851: }

855: PetscErrorCode VecView_Seq_Draw(Vec xin,PetscViewer v)
856: {
857:   PetscErrorCode    ierr;
858:   PetscDraw         draw;
859:   PetscBool         isnull;
860:   PetscViewerFormat format;

863:   PetscViewerDrawGetDraw(v,0,&draw);
864:   PetscDrawIsNull(draw,&isnull); if (isnull) return(0);

866:   PetscViewerGetFormat(v,&format);
867:   /*
868:      Currently it only supports drawing to a line graph */
869:   if (format != PETSC_VIEWER_DRAW_LG) {
870:     PetscViewerPushFormat(v,PETSC_VIEWER_DRAW_LG);
871:   }
872:   VecView_Seq_Draw_LG(xin,v);
873:   if (format != PETSC_VIEWER_DRAW_LG) {
874:     PetscViewerPopFormat(v);
875:   }
876:   return(0);
877: }

881: PetscErrorCode VecView_Seq_Binary(Vec xin,PetscViewer viewer)
882: {
883:   PetscErrorCode    ierr;
884:   int               fdes;
885:   PetscInt          n = xin->map->n,classid=VEC_FILE_CLASSID;
886:   FILE              *file;
887:   const PetscScalar *xv;
888: #if defined(PETSC_HAVE_MPIIO)
889:   PetscBool         isMPIIO;
890: #endif
891:   PetscBool         skipHeader;
892:   PetscViewerFormat format;

895:   /* Write vector header */
896:   PetscViewerBinaryGetSkipHeader(viewer,&skipHeader);
897:   if (!skipHeader) {
898:     PetscViewerBinaryWrite(viewer,&classid,1,PETSC_INT,PETSC_FALSE);
899:     PetscViewerBinaryWrite(viewer,&n,1,PETSC_INT,PETSC_FALSE);
900:   }

902:   /* Write vector contents */
903: #if defined(PETSC_HAVE_MPIIO)
904:   PetscViewerBinaryGetMPIIO(viewer,&isMPIIO);
905:   if (!isMPIIO) {
906: #endif
907:     PetscViewerBinaryGetDescriptor(viewer,&fdes);
908:     VecGetArrayRead(xin,&xv);
909:     PetscBinaryWrite(fdes,(void*)xv,n,PETSC_SCALAR,PETSC_FALSE);
910:     VecRestoreArrayRead(xin,&xv);
911:     PetscViewerGetFormat(viewer,&format);
912:     if (format == PETSC_VIEWER_BINARY_MATLAB) {
913:       MPI_Comm   comm;
914:       FILE       *info;
915:       const char *name;

917:       PetscObjectGetName((PetscObject)xin,&name);
918:       PetscObjectGetComm((PetscObject)viewer,&comm);
919:       PetscViewerBinaryGetInfoPointer(viewer,&info);
920:       PetscFPrintf(comm,info,"%%--- begin code written by PetscViewerBinary for MATLAB format ---%\n");
921:       PetscFPrintf(comm,info,"%%$$ Set.%s = PetscBinaryRead(fd);\n",name);
922:       PetscFPrintf(comm,info,"%%--- end code written by PetscViewerBinary for MATLAB format ---%\n\n");
923:     }
924: #if defined(PETSC_HAVE_MPIIO)
925:   } else {
926:     MPI_Offset   off;
927:     MPI_File     mfdes;
928:     PetscMPIInt  lsize;

930:     PetscMPIIntCast(n,&lsize);
931:     PetscViewerBinaryGetMPIIODescriptor(viewer,&mfdes);
932:     PetscViewerBinaryGetMPIIOOffset(viewer,&off);
933:     MPI_File_set_view(mfdes,off,MPIU_SCALAR,MPIU_SCALAR,(char*)"native",MPI_INFO_NULL);
934:     VecGetArrayRead(xin,&xv);
935:     MPIU_File_write_all(mfdes,(void*)xv,lsize,MPIU_SCALAR,MPI_STATUS_IGNORE);
936:     VecRestoreArrayRead(xin,&xv);
937:     PetscViewerBinaryAddMPIIOOffset(viewer,n*sizeof(PetscScalar));
938:   }
939: #endif

941:   PetscViewerBinaryGetInfoPointer(viewer,&file);
942:   if (file) {
943:     if (((PetscObject)xin)->prefix) {
944:       PetscFPrintf(PETSC_COMM_SELF,file,"-%s_vecload_block_size %D\n",((PetscObject)xin)->prefix,xin->map->bs);
945:     } else {
946:       PetscFPrintf(PETSC_COMM_SELF,file,"-vecload_block_size %D\n",xin->map->bs);
947:     }
948:   }
949:   return(0);
950: }

952: #if defined(PETSC_HAVE_MATLAB_ENGINE)
953: #include <petscmatlab.h>
954: #include <mat.h>   /* MATLAB include file */
957: PetscErrorCode VecView_Seq_Matlab(Vec vec,PetscViewer viewer)
958: {
959:   PetscErrorCode    ierr;
960:   PetscInt          n;
961:   const PetscScalar *array;

964:   VecGetLocalSize(vec,&n);
965:   PetscObjectName((PetscObject)vec);
966:   VecGetArrayRead(vec,&array);
967:   PetscViewerMatlabPutArray(viewer,n,1,array,((PetscObject)vec)->name);
968:   VecRestoreArrayRead(vec,&array);
969:   return(0);
970: }
971: #endif

975: PetscErrorCode VecView_Seq(Vec xin,PetscViewer viewer)
976: {
978:   PetscBool      isdraw,iascii,issocket,isbinary;
979: #if defined(PETSC_HAVE_MATHEMATICA)
980:   PetscBool      ismathematica;
981: #endif
982: #if defined(PETSC_HAVE_MATLAB_ENGINE)
983:   PetscBool      ismatlab;
984: #endif
985: #if defined(PETSC_HAVE_HDF5)
986:   PetscBool      ishdf5;
987: #endif

990:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERDRAW,&isdraw);
991:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);
992:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERSOCKET,&issocket);
993:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERBINARY,&isbinary);
994: #if defined(PETSC_HAVE_MATHEMATICA)
995:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERMATHEMATICA,&ismathematica);
996: #endif
997: #if defined(PETSC_HAVE_HDF5)
998:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERHDF5,&ishdf5);
999: #endif
1000: #if defined(PETSC_HAVE_MATLAB_ENGINE)
1001:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERMATLAB,&ismatlab);
1002: #endif

1004:   if (isdraw) {
1005:     VecView_Seq_Draw(xin,viewer);
1006:   } else if (iascii) {
1007:     VecView_Seq_ASCII(xin,viewer);
1008:   } else if (isbinary) {
1009:     VecView_Seq_Binary(xin,viewer);
1010: #if defined(PETSC_HAVE_MATHEMATICA)
1011:   } else if (ismathematica) {
1012:     PetscViewerMathematicaPutVector(viewer,xin);
1013: #endif
1014: #if defined(PETSC_HAVE_HDF5)
1015:   } else if (ishdf5) {
1016:     VecView_MPI_HDF5(xin,viewer); /* Reusing VecView_MPI_HDF5 ... don't want code duplication*/
1017: #endif
1018: #if defined(PETSC_HAVE_MATLAB_ENGINE)
1019:   } else if (ismatlab) {
1020:     VecView_Seq_Matlab(xin,viewer);
1021: #endif
1022:   }
1023:   return(0);
1024: }

1028: PetscErrorCode VecGetValues_Seq(Vec xin,PetscInt ni,const PetscInt ix[],PetscScalar y[])
1029: {
1030:   const PetscScalar *xx;
1031:   PetscInt          i;
1032:   PetscErrorCode    ierr;

1035:   VecGetArrayRead(xin,&xx);
1036:   for (i=0; i<ni; i++) {
1037:     if (xin->stash.ignorenegidx && ix[i] < 0) continue;
1038: #if defined(PETSC_USE_DEBUG)
1039:     if (ix[i] < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Out of range index value %D cannot be negative",ix[i]);
1040:     if (ix[i] >= xin->map->n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Out of range index value %D to large maximum allowed %D",ix[i],xin->map->n);
1041: #endif
1042:     y[i] = xx[ix[i]];
1043:   }
1044:   VecRestoreArrayRead(xin,&xx);
1045:   return(0);
1046: }

1050: PetscErrorCode VecSetValues_Seq(Vec xin,PetscInt ni,const PetscInt ix[],const PetscScalar y[],InsertMode m)
1051: {
1052:   PetscScalar    *xx;
1053:   PetscInt       i;

1057:   VecGetArray(xin,&xx);
1058:   if (m == INSERT_VALUES) {
1059:     for (i=0; i<ni; i++) {
1060:       if (xin->stash.ignorenegidx && ix[i] < 0) continue;
1061: #if defined(PETSC_USE_DEBUG)
1062:       if (ix[i] < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Out of range index value %D cannot be negative",ix[i]);
1063:       if (ix[i] >= xin->map->n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Out of range index value %D maximum %D",ix[i],xin->map->n);
1064: #endif
1065:       xx[ix[i]] = y[i];
1066:     }
1067:   } else {
1068:     for (i=0; i<ni; i++) {
1069:       if (xin->stash.ignorenegidx && ix[i] < 0) continue;
1070: #if defined(PETSC_USE_DEBUG)
1071:       if (ix[i] < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Out of range index value %D cannot be negative",ix[i]);
1072:       if (ix[i] >= xin->map->n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Out of range index value %D maximum %D",ix[i],xin->map->n);
1073: #endif
1074:       xx[ix[i]] += y[i];
1075:     }
1076:   }
1077:   VecRestoreArray(xin,&xx);
1078:   return(0);
1079: }

1083: PetscErrorCode VecSetValuesBlocked_Seq(Vec xin,PetscInt ni,const PetscInt ix[],const PetscScalar yin[],InsertMode m)
1084: {
1085:   PetscScalar    *xx,*y = (PetscScalar*)yin;
1086:   PetscInt       i,bs = xin->map->bs,start,j;

1089:   /*
1090:        For optimization could treat bs = 2, 3, 4, 5 as special cases with loop unrolling
1091:   */
1093:   VecGetArray(xin,&xx);
1094:   if (m == INSERT_VALUES) {
1095:     for (i=0; i<ni; i++) {
1096:       start = bs*ix[i];
1097:       if (start < 0) continue;
1098: #if defined(PETSC_USE_DEBUG)
1099:       if (start >= xin->map->n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Out of range index value %D maximum %D",start,xin->map->n);
1100: #endif
1101:       for (j=0; j<bs; j++) xx[start+j] = y[j];
1102:       y += bs;
1103:     }
1104:   } else {
1105:     for (i=0; i<ni; i++) {
1106:       start = bs*ix[i];
1107:       if (start < 0) continue;
1108: #if defined(PETSC_USE_DEBUG)
1109:       if (start >= xin->map->n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Out of range index value %D maximum %D",start,xin->map->n);
1110: #endif
1111:       for (j=0; j<bs; j++) xx[start+j] += y[j];
1112:       y += bs;
1113:     }
1114:   }
1115:   VecRestoreArray(xin,&xx);
1116:   return(0);
1117: }

1121: PetscErrorCode VecDestroy_Seq(Vec v)
1122: {
1123:   Vec_Seq        *vs = (Vec_Seq*)v->data;

1127: #if defined(PETSC_USE_LOG)
1128:   PetscLogObjectState((PetscObject)v,"Length=%D",v->map->n);
1129: #endif
1130:   PetscFree(vs->array_allocated);
1131:   PetscFree(v->data);
1132:   return(0);
1133: }

1137: PetscErrorCode VecSetOption_Seq(Vec v,VecOption op,PetscBool flag)
1138: {
1140:   if (op == VEC_IGNORE_NEGATIVE_INDICES) v->stash.ignorenegidx = flag;
1141:   return(0);
1142: }

1146: PetscErrorCode VecDuplicate_Seq(Vec win,Vec *V)
1147: {

1151:   VecCreate(PetscObjectComm((PetscObject)win),V);
1152:   PetscObjectSetPrecision((PetscObject)*V,((PetscObject)win)->precision);
1153:   VecSetSizes(*V,win->map->n,win->map->n);
1154:   VecSetType(*V,((PetscObject)win)->type_name);
1155:   PetscLayoutReference(win->map,&(*V)->map);
1156:   PetscObjectListDuplicate(((PetscObject)win)->olist,&((PetscObject)(*V))->olist);
1157:   PetscFunctionListDuplicate(((PetscObject)win)->qlist,&((PetscObject)(*V))->qlist);

1159:   (*V)->ops->view          = win->ops->view;
1160:   (*V)->stash.ignorenegidx = win->stash.ignorenegidx;
1161:   return(0);
1162: }

1164: static struct _VecOps DvOps = {VecDuplicate_Seq, /* 1 */
1165:                                VecDuplicateVecs_Default,
1166:                                VecDestroyVecs_Default,
1167:                                VecDot_Seq,
1168:                                VecMDot_Seq,
1169:                                VecNorm_Seq,
1170:                                VecTDot_Seq,
1171:                                VecMTDot_Seq,
1172:                                VecScale_Seq,
1173:                                VecCopy_Seq, /* 10 */
1174:                                VecSet_Seq,
1175:                                VecSwap_Seq,
1176:                                VecAXPY_Seq,
1177:                                VecAXPBY_Seq,
1178:                                VecMAXPY_Seq,
1179:                                VecAYPX_Seq,
1180:                                VecWAXPY_Seq,
1181:                                VecAXPBYPCZ_Seq,
1182:                                VecPointwiseMult_Seq,
1183:                                VecPointwiseDivide_Seq,
1184:                                VecSetValues_Seq, /* 20 */
1185:                                0,0,
1186:                                0,
1187:                                VecGetSize_Seq,
1188:                                VecGetSize_Seq,
1189:                                0,
1190:                                VecMax_Seq,
1191:                                VecMin_Seq,
1192:                                VecSetRandom_Seq,
1193:                                VecSetOption_Seq, /* 30 */
1194:                                VecSetValuesBlocked_Seq,
1195:                                VecDestroy_Seq,
1196:                                VecView_Seq,
1197:                                VecPlaceArray_Seq,
1198:                                VecReplaceArray_Seq,
1199:                                VecDot_Seq,
1200:                                VecTDot_Seq,
1201:                                VecNorm_Seq,
1202:                                VecMDot_Seq,
1203:                                VecMTDot_Seq, /* 40 */
1204:                                VecLoad_Default,
1205:                                VecReciprocal_Default,
1206:                                VecConjugate_Seq,
1207:                                0,
1208:                                0,
1209:                                VecResetArray_Seq,
1210:                                0,
1211:                                VecMaxPointwiseDivide_Seq,
1212:                                VecPointwiseMax_Seq,
1213:                                VecPointwiseMaxAbs_Seq,
1214:                                VecPointwiseMin_Seq,
1215:                                VecGetValues_Seq,
1216:                                0,
1217:                                0,
1218:                                0,
1219:                                0,
1220:                                0,
1221:                                0,
1222:                                VecStrideGather_Default,
1223:                                VecStrideScatter_Default,
1224:                                0,
1225:                                0,
1226:                                0,
1227:                                0,
1228:                                0
1229: };


1232: /*
1233:       This is called by VecCreate_Seq() (i.e. VecCreateSeq()) and VecCreateSeqWithArray()
1234: */
1237: PetscErrorCode VecCreate_Seq_Private(Vec v,const PetscScalar array[])
1238: {
1239:   Vec_Seq        *s;

1243:   PetscNewLog(v,Vec_Seq,&s);
1244:   PetscMemcpy(v->ops,&DvOps,sizeof(DvOps));

1246:   v->data            = (void*)s;
1247:   v->petscnative     = PETSC_TRUE;
1248:   s->array           = (PetscScalar*)array;
1249:   s->array_allocated = 0;

1251:   PetscLayoutSetUp(v->map);
1252:   PetscObjectChangeTypeName((PetscObject)v,VECSEQ);
1253: #if defined(PETSC_HAVE_MATLAB_ENGINE)
1254:   PetscObjectComposeFunction((PetscObject)v,"PetscMatlabEnginePut_C",VecMatlabEnginePut_Default);
1255:   PetscObjectComposeFunction((PetscObject)v,"PetscMatlabEngineGet_C",VecMatlabEngineGet_Default);
1256: #endif
1257: #if defined(PETSC_USE_MIXED_PRECISION)
1258:   ((PetscObject)v)->precision = (PetscPrecision)sizeof(PetscReal);
1259: #endif
1260:   return(0);
1261: }

1265: /*@C
1266:    VecCreateSeqWithArray - Creates a standard,sequential array-style vector,
1267:    where the user provides the array space to store the vector values.

1269:    Collective on MPI_Comm

1271:    Input Parameter:
1272: +  comm - the communicator, should be PETSC_COMM_SELF
1273: .  n - the vector length
1274: -  array - memory where the vector elements are to be stored.

1276:    Output Parameter:
1277: .  V - the vector

1279:    Notes:
1280:    Use VecDuplicate() or VecDuplicateVecs() to form additional vectors of the
1281:    same type as an existing vector.

1283:    If the user-provided array is NULL, then VecPlaceArray() can be used
1284:    at a later stage to SET the array for storing the vector values.

1286:    PETSc does NOT free the array when the vector is destroyed via VecDestroy().
1287:    The user should not free the array until the vector is destroyed.

1289:    Level: intermediate

1291:    Concepts: vectors^creating with array

1293: .seealso: VecCreateMPIWithArray(), VecCreate(), VecDuplicate(), VecDuplicateVecs(),
1294:           VecCreateGhost(), VecCreateSeq(), VecPlaceArray()
1295: @*/
1296: PetscErrorCode  VecCreateSeqWithArray(MPI_Comm comm,PetscInt bs,PetscInt n,const PetscScalar array[],Vec *V)
1297: {
1299:   PetscMPIInt    size;

1302:   VecCreate(comm,V);
1303:   VecSetSizes(*V,n,n);
1304:   VecSetBlockSize(*V,bs);
1305:   MPI_Comm_size(comm,&size);
1306:   if (size > 1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Cannot create VECSEQ on more than one process");
1307:   VecCreate_Seq_Private(*V,array);
1308:   return(0);
1309: }