xref: /petsc/src/sys/tests/ex10.c (revision c4762a1b19cd2af06abeed90e8f9d34fb975dd94)
1*c4762a1bSJed Brown 
2*c4762a1bSJed Brown static char help[] = "Tests PetscArraymove()/PetscMemmove()\n";
3*c4762a1bSJed Brown 
4*c4762a1bSJed Brown #include <petscsys.h>
5*c4762a1bSJed Brown 
6*c4762a1bSJed Brown int main(int argc,char **argv)
7*c4762a1bSJed Brown {
8*c4762a1bSJed Brown   PetscInt       i,*a,*b;
9*c4762a1bSJed Brown   PetscErrorCode ierr;
10*c4762a1bSJed Brown 
11*c4762a1bSJed Brown   ierr = PetscInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr;
12*c4762a1bSJed Brown 
13*c4762a1bSJed Brown   ierr = PetscMalloc1(10,&a);CHKERRQ(ierr);
14*c4762a1bSJed Brown   ierr = PetscMalloc1(20,&b);CHKERRQ(ierr);
15*c4762a1bSJed Brown 
16*c4762a1bSJed Brown   /*
17*c4762a1bSJed Brown       Nonoverlapping regions
18*c4762a1bSJed Brown   */
19*c4762a1bSJed Brown   for (i=0; i<20; i++) b[i] = i;
20*c4762a1bSJed Brown   ierr = PetscArraymove(a,b,10);CHKERRQ(ierr);
21*c4762a1bSJed Brown   ierr = PetscIntView(10,a,NULL);CHKERRQ(ierr);
22*c4762a1bSJed Brown 
23*c4762a1bSJed Brown   ierr = PetscFree(a);CHKERRQ(ierr);
24*c4762a1bSJed Brown 
25*c4762a1bSJed Brown   /*
26*c4762a1bSJed Brown      |        |                |       |
27*c4762a1bSJed Brown      b        a               b+15    b+20
28*c4762a1bSJed Brown                               a+10    a+15
29*c4762a1bSJed Brown   */
30*c4762a1bSJed Brown   a    = b + 5;
31*c4762a1bSJed Brown   ierr = PetscArraymove(a,b,15);CHKERRQ(ierr);
32*c4762a1bSJed Brown   ierr = PetscIntView(15,a,NULL);CHKERRQ(ierr);
33*c4762a1bSJed Brown   ierr = PetscFree(b);CHKERRQ(ierr);
34*c4762a1bSJed Brown 
35*c4762a1bSJed Brown   /*
36*c4762a1bSJed Brown      |       |                    |       |
37*c4762a1bSJed Brown      a       b                   a+20   a+25
38*c4762a1bSJed Brown                                         b+20
39*c4762a1bSJed Brown   */
40*c4762a1bSJed Brown   ierr = PetscMalloc1(25,&a);CHKERRQ(ierr);
41*c4762a1bSJed Brown   b    = a + 5;
42*c4762a1bSJed Brown   for (i=0; i<20; i++) b[i] = i;
43*c4762a1bSJed Brown   ierr = PetscArraymove(a,b,20);CHKERRQ(ierr);
44*c4762a1bSJed Brown   ierr = PetscIntView(20,a,NULL);CHKERRQ(ierr);
45*c4762a1bSJed Brown   ierr = PetscFree(a);CHKERRQ(ierr);
46*c4762a1bSJed Brown 
47*c4762a1bSJed Brown   ierr = PetscFinalize();
48*c4762a1bSJed Brown   return ierr;
49*c4762a1bSJed Brown }
50*c4762a1bSJed Brown 
51*c4762a1bSJed Brown 
52*c4762a1bSJed Brown 
53*c4762a1bSJed Brown /*TEST
54*c4762a1bSJed Brown 
55*c4762a1bSJed Brown    test:
56*c4762a1bSJed Brown 
57*c4762a1bSJed Brown TEST*/
58