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