[petsc-users] AddressSanitizer: attempting free on address which was not malloc()-ed

Jed Brown jed at jedbrown.org
Sun Mar 3 16:33:04 CST 2019


The compiler doesn't know anything special about PetscFinalize.
Destructors are called after all executable statements in their scope.
So you need the extra scoping if the destructor should be called
earlier.

Yuyun Yang <yyang85 at stanford.edu> writes:

> Oh interesting, so I need to add those extra brackets around my class object and function calls. I thought the destructor is automatically at Finalize.
>
> Thanks!
> Yuyun
>
> -----Original Message-----
> From: Jed Brown <jed at jedbrown.org> 
> Sent: Sunday, March 3, 2019 2:19 PM
> To: Yuyun Yang <yyang85 at stanford.edu>; Matthew Knepley <knepley at gmail.com>
> Cc: petsc-users at mcs.anl.gov
> Subject: Re: [petsc-users] AddressSanitizer: attempting free on address which was not malloc()-ed
>
> If you run this with MPICH, it prints
>
>   Attempting to use an MPI routine after finalizing MPICH
>
> You need to ensure that the C++ class destructor is called before PetscFinalize.  For example, like this:
>
> diff --git i/test_domain.cpp w/test_domain.cpp index 0cfe22f..23545f2 100644
> --- i/test_domain.cpp
> +++ w/test_domain.cpp
> @@ -8,11 +8,12 @@ int main(int argc, char **argv) {
>    PetscErrorCode ierr = 0;
>    PetscInitialize(&argc, &argv, NULL, NULL);
>  
> -  Domain d;
> +  {
> +    Domain d;
>  
> -  ierr = d.setFields(); CHKERRQ(ierr);
> -  ierr = d.setScatters(); CHKERRQ(ierr);
> -  
> +    ierr = d.setFields(); CHKERRQ(ierr);
> +    ierr = d.setScatters(); CHKERRQ(ierr);  }
>    PetscFinalize();
>    return ierr;
>  }
>
>
> Yuyun Yang via petsc-users <petsc-users at mcs.anl.gov> writes:
>
>> Yes, please see the attached files for a minimal example. Thanks a lot!
>>
>> Best regards,
>> Yuyun
>>
>> From: Matthew Knepley <knepley at gmail.com>
>> Sent: Sunday, March 3, 2019 12:46 PM
>> To: Yuyun Yang <yyang85 at stanford.edu>
>> Cc: Zhang, Junchao <jczhang at mcs.anl.gov>; petsc-users at mcs.anl.gov
>> Subject: Re: [petsc-users] AddressSanitizer: attempting free on 
>> address which was not malloc()-ed
>>
>> On Sun, Mar 3, 2019 at 3:05 PM Yuyun Yang <yyang85 at stanford.edu<mailto:yyang85 at stanford.edu>> wrote:
>> Actually, I tried just creating a domain object (since the address sanitizer was complaining about that code to start with). Simply creating that object gave me a core dump, so I suppose the issue must be there. I got the following message when running the code with -objects_dump flag on the command line, but couldn’t find a problem with the code (I’ve attached it here with only the relevant functions).
>>
>> I think what we are going to need from you is a minimal example that 
>> has the error. I am guessing you have a logic bug in the C++, which we cannot debug by looking.
>>
>>   Thanks,
>>
>>      Matt
>>
>> Thanks a lot for your help!
>>
>> The following objects were never freed
>> -----------------------------------------
>> [0] Vec seq y
>>       [0]  VecCreate() in 
>> /home/yyy910805/petsc/src/vec/vec/interface/veccreate.c
>> [0] Vec seq Vec_0x84000000_0
>>       [0]  VecCreate() in 
>> /home/yyy910805/petsc/src/vec/vec/interface/veccreate.c
>> [0] Vec seq Vec_0x84000000_1
>>       [0]  VecCreate() in 
>> /home/yyy910805/petsc/src/vec/vec/interface/veccreate.c
>> [0] VecScatter seq VecScatter_0x84000000_2
>>       [0]  VecScatterCreate() in 
>> /home/yyy910805/petsc/src/vec/vscat/interface/vscreate.c
>> [0] VecScatter seq VecScatter_0x84000000_3
>>       [0]  VecScatterCreate() in 
>> /home/yyy910805/petsc/src/vec/vscat/interface/vscreate.c
>> [0] VecScatter seq VecScatter_0x84000000_4
>>       [0]  VecScatterCreate() in 
>> /home/yyy910805/petsc/src/vec/vscat/interface/vscreate.c
>> [0] VecScatter seq VecScatter_0x84000000_5
>>       [0]  VecScatterCreate() in 
>> /home/yyy910805/petsc/src/vec/vscat/interface/vscreate.c
>> Attempting to use an MPI routine after finalizing MPICH
>>
>> ----------------------------------------------------------------------
>> ----------------------------------------------------------------------
>> -----
>> From: Matthew Knepley <knepley at gmail.com<mailto:knepley at gmail.com>>
>> Sent: Sunday, March 3, 2019 11:28 AM
>> To: Yuyun Yang <yyang85 at stanford.edu<mailto:yyang85 at stanford.edu>>
>> Cc: Zhang, Junchao <jczhang at mcs.anl.gov<mailto:jczhang at mcs.anl.gov>>; 
>> petsc-users at mcs.anl.gov<mailto:petsc-users at mcs.anl.gov>
>> Subject: Re: [petsc-users] AddressSanitizer: attempting free on 
>> address which was not malloc()-ed
>>
>> On Sun, Mar 3, 2019 at 1:03 PM Yuyun Yang via petsc-users <petsc-users at mcs.anl.gov<mailto:petsc-users at mcs.anl.gov>> wrote:
>> I tried compiling without the sanitizer and running on valgrind. Got a bunch of errors “Uninitialised value was created by a stack allocation at 0x41B280: ComputeVel_qd::computeVel(double*, double, int&, int)”.
>>
>> There is no memory management code here, so other parts of the code must be relevant.
>>
>>   Thanks,
>>
>>     Matt
>>
>>
>> HEAP SUMMARY:
>> ==74==     in use at exit: 96,637 bytes in 91 blocks
>> ==74==   total heap usage: 47,774 allocs, 47,522 frees, 308,253,653 bytes allocated
>> LEAK SUMMARY:
>> ==74==    definitely lost: 0 bytes in 0 blocks
>> ==74==    indirectly lost: 0 bytes in 0 blocks
>> ==74==      possibly lost: 0 bytes in 0 blocks
>> ==74==    still reachable: 96,637 bytes in 91 blocks
>> ==74==         suppressed: 0 bytes in 0 blocks
>>
>> The error is located in the attached code (I’ve extracted only the relevant functions), but I couldn’t figure out what is wrong. Is this causing the memory corruption/double free error that happens when I execute the code?
>>
>> Thanks a lot for your help.
>>
>> Best regards,
>> Yuyun
>>
>> From: Zhang, Junchao <jczhang at mcs.anl.gov<mailto:jczhang at mcs.anl.gov>>
>> Sent: Friday, March 1, 2019 7:36 AM
>> To: Yuyun Yang <yyang85 at stanford.edu<mailto:yyang85 at stanford.edu>>
>> Subject: Re: [petsc-users] AddressSanitizer: attempting free on 
>> address which was not malloc()-ed
>>
>>
>> On Fri, Mar 1, 2019 at 1:02 AM Yuyun Yang <yyang85 at stanford.edu<mailto:yyang85 at stanford.edu>> wrote:
>> Actually, I also saw a line at the beginning of valgrind saying "shadow memory range interleaves with an existing memory mapping. ASan cannot proceed properly. ABORTING." I guess the code didn't really run through valgrind since it aborted. Should I remove the address sanitizer flag when compiling?
>> From the message, it seems ASan (not valgrind) aborted. You can try to compile without sanitizer and then run with valgrind. If no problem, then it is probably a sanitizer issue.
>>
>>
>> Get Outlook for iOS<https://aka.ms/o0ukef> 
>> ________________________________
>> From: Yuyun Yang
>> Sent: Thursday, February 28, 2019 10:54:57 PM
>> To: Zhang, Junchao
>> Subject: Re: [petsc-users] AddressSanitizer: attempting free on 
>> address which was not malloc()-ed
>>
>> Hmm, still getting the same error from address sanitizer even though valgrind shows no errors and no leaks are possible.
>>
>> Should I ignore that error? My results did run alright.
>>
>> Best,
>> Yuyun
>>
>> Get Outlook for iOS<https://aka.ms/o0ukef> 
>> ________________________________
>> From: Zhang, Junchao <jczhang at mcs.anl.gov<mailto:jczhang at mcs.anl.gov>>
>> Sent: Wednesday, February 27, 2019 8:27:17 PM
>> To: Yuyun Yang
>> Cc: Matthew Knepley; 
>> petsc-users at mcs.anl.gov<mailto:petsc-users at mcs.anl.gov>
>> Subject: Re: [petsc-users] AddressSanitizer: attempting free on 
>> address which was not malloc()-ed
>>
>> Try the following to see if you can catch the bug easily: 1) Get error 
>> code for each petsc function and check it with CHKERRQ; 2) Link your 
>> code with a petsc library with debugging enabled (configured with 
>> --with-debugging=1); 3) Run your code with valgrind
>>
>> --Junchao Zhang
>>
>>
>> On Wed, Feb 27, 2019 at 9:04 PM Yuyun Yang <yyang85 at stanford.edu<mailto:yyang85 at stanford.edu>> wrote:
>> Hi Junchao,
>>
>> This code actually involves a lot of classes and is pretty big. Might be an overkill for me to send everything to you. I'd like to know if I see this sort of error message, which points to this domain file, is it possible that the problem happens in another file (whose operations are linked to this one)? If so, I'll debug a little more and maybe send you more useful information later.
>>
>> Best regards,
>> Yuyun
>>
>> Get Outlook for iOS<https://aka.ms/o0ukef> 
>> ________________________________
>> From: Zhang, Junchao <jczhang at mcs.anl.gov<mailto:jczhang at mcs.anl.gov>>
>> Sent: Wednesday, February 27, 2019 6:24:13 PM
>> To: Yuyun Yang
>> Cc: Matthew Knepley; 
>> petsc-users at mcs.anl.gov<mailto:petsc-users at mcs.anl.gov>
>> Subject: Re: [petsc-users] AddressSanitizer: attempting free on 
>> address which was not malloc()-ed
>>
>> Could you provide a compilable and runnable test so I can try it?
>> --Junchao Zhang
>>
>>
>> On Wed, Feb 27, 2019 at 7:34 PM Yuyun Yang <yyang85 at stanford.edu<mailto:yyang85 at stanford.edu>> wrote:
>> Thanks, I fixed that, but I’m not actually calling the testScatters() function in my implementation (in the constructor, the only functions I called are setFields and setScatters). So the problem couldn’t have been that?
>>
>> Best,
>> Yuyun
>>
>> From: Zhang, Junchao <jczhang at mcs.anl.gov<mailto:jczhang at mcs.anl.gov>>
>> Sent: Wednesday, February 27, 2019 10:50 AM
>> To: Yuyun Yang <yyang85 at stanford.edu<mailto:yyang85 at stanford.edu>>
>> Cc: Matthew Knepley <knepley at gmail.com<mailto:knepley at gmail.com>>; 
>> petsc-users at mcs.anl.gov<mailto:petsc-users at mcs.anl.gov>
>> Subject: Re: [petsc-users] AddressSanitizer: attempting free on 
>> address which was not malloc()-ed
>>
>>
>> On Wed, Feb 27, 2019 at 10:41 AM Yuyun Yang via petsc-users <petsc-users at mcs.anl.gov<mailto:petsc-users at mcs.anl.gov>> wrote:
>> I called VecDestroy() in the destructor for this object – is that not the right way to do it?
>> In Domain::testScatters(), you have many VecDuplicate(,&out), You need 
>> to VecDestroy(&out) before doing new VecDuplicate(,&out); How do I implement CHECK ALL RETURN CODES?
>> For each PETSc function, do ierr = ...;  CHKERRQ(ierr);
>>
>> From: Matthew Knepley <knepley at gmail.com<mailto:knepley at gmail.com>>
>> Sent: Wednesday, February 27, 2019 7:24 AM
>> To: Yuyun Yang <yyang85 at stanford.edu<mailto:yyang85 at stanford.edu>>
>> Cc: petsc-users at mcs.anl.gov<mailto:petsc-users at mcs.anl.gov>
>> Subject: Re: [petsc-users] AddressSanitizer: attempting free on 
>> address which was not malloc()-ed
>>
>> You call VecDuplicate() a bunch, but VecDestroy() only once in the bottom function. This is wrong.
>> Also, CHECK ALL RETURN CODES. This is the fastest way to find errors.
>>
>>    Matt
>>
>> On Wed, Feb 27, 2019 at 2:06 AM Yuyun Yang via petsc-users <petsc-users at mcs.anl.gov<mailto:petsc-users at mcs.anl.gov>> wrote:
>> Hello team,
>>
>> I ran into the address sanitizer error that I hope you could help me with. I don’t really know what’s wrong with the way the code frees memory. The relevant code file is attached. The line number following domain.cpp specifically referenced to the vector _q, which seems a little odd, since some other vectors are constructed and freed the same way.
>>
>> ==1719==ERROR: AddressSanitizer: attempting free on address which was not malloc()-ed: 0x61f0000076c0 in thread T0
>>     #0 0x7fbf195282ca in __interceptor_free (/usr/lib/x86_64-linux-gnu/libasan.so.2+0x982ca)
>>     #1 0x7fbf1706f895 in PetscFreeAlign /home/yyy910805/petsc/src/sys/memory/mal.c:87
>>     #2 0x7fbf1731a898 in VecDestroy_Seq /home/yyy910805/petsc/src/vec/vec/impls/seq/bvec2.c:788
>>     #3 0x7fbf1735f795 in VecDestroy /home/yyy910805/petsc/src/vec/vec/interface/vector.c:408
>>     #4 0x40dd0a in Domain::~Domain() /home/yyy910805/scycle/source/domain.cpp:132
>>     #5 0x40b479 in main /home/yyy910805/scycle/source/main.cpp:242
>>     #6 0x7fbf14d2082f in __libc_start_main (/lib/x86_64-linux-gnu/libc.so.6+0x2082f)
>>     #7 0x4075d8 in _start 
>> (/home/yyy910805/scycle/source/main+0x4075d8)
>>
>> 0x61f0000076c0 is located 1600 bytes inside of 3220-byte region 
>> [0x61f000007080,0x61f000007d14) allocated by thread T0 here:
>>     #0 0x7fbf19528b32 in __interceptor_memalign (/usr/lib/x86_64-linux-gnu/libasan.so.2+0x98b32)
>>     #1 0x7fbf1706f7e0 in PetscMallocAlign /home/yyy910805/petsc/src/sys/memory/mal.c:41
>>     #2 0x7fbf17073022 in PetscTrMallocDefault /home/yyy910805/petsc/src/sys/memory/mtr.c:183
>>     #3 0x7fbf170710a1 in PetscMallocA /home/yyy910805/petsc/src/sys/memory/mal.c:397
>>     #4 0x7fbf17326fb0 in VecCreate_Seq /home/yyy910805/petsc/src/vec/vec/impls/seq/bvec3.c:35
>>     #5 0x7fbf1736f560 in VecSetType /home/yyy910805/petsc/src/vec/vec/interface/vecreg.c:51
>>     #6 0x7fbf1731afae in VecDuplicate_Seq /home/yyy910805/petsc/src/vec/vec/impls/seq/bvec2.c:807
>>     #7 0x7fbf1735eff7 in VecDuplicate /home/yyy910805/petsc/src/vec/vec/interface/vector.c:379
>>     #8 0x4130de in Domain::setFields() /home/yyy910805/scycle/source/domain.cpp:431
>>     #9 0x40c60a in Domain::Domain(char const*) /home/yyy910805/scycle/source/domain.cpp:57
>>     #10 0x40b433 in main /home/yyy910805/scycle/source/main.cpp:242
>>     #11 0x7fbf14d2082f in __libc_start_main 
>> (/lib/x86_64-linux-gnu/libc.so.6+0x2082f)
>>
>> SUMMARY: AddressSanitizer: bad-free ??:0 __interceptor_free 
>> ==1719==ABORTING
>>
>> Thanks very much!
>> Yuyun
>>
>>
>> --
>> What most experimenters take for granted before they begin their experiments is infinitely more interesting than any results to which their experiments lead.
>> -- Norbert Wiener
>>
>> https://www.cse.buffalo.edu/~knepley/<http://www.cse.buffalo.edu/~knep
>> ley/>
>>
>>
>> --
>> What most experimenters take for granted before they begin their experiments is infinitely more interesting than any results to which their experiments lead.
>> -- Norbert Wiener
>>
>> https://www.cse.buffalo.edu/~knepley/<http://www.cse.buffalo.edu/~knep
>> ley/>
>>
>>
>> --
>> What most experimenters take for granted before they begin their experiments is infinitely more interesting than any results to which their experiments lead.
>> -- Norbert Wiener
>>
>> https://www.cse.buffalo.edu/~knepley/<http://www.cse.buffalo.edu/~knep
>> ley/>
>> #include "domain.hpp"
>>
>> using namespace std;
>>
>> Domain::Domain()
>>   :  _sbpType("mfc_coordTrans"),_Ny(201),_Nz(1),_Ly(30),_Lz(30),
>>      _q(NULL),_r(NULL),_y(NULL),_z(NULL),_y0(NULL),_z0(NULL),_dq(-1),_dr(-1),
>>      _bCoordTrans(5)
>> {
>>   if (_Ny > 1) {
>>     _dq = 1.0 / (_Ny - 1.0);
>>   }
>>   else {
>>     _dq = 1;
>>   }
>>
>>   if (_Nz > 1) {
>>     _dr = 1.0 / (_Nz - 1.0);
>>   }
>>   else {
>>     _dr = 1;
>>   }
>> }
>>
>> // destructor
>> Domain::~Domain()
>> {
>>   // free memory
>>   VecDestroy(&_q);
>>   VecDestroy(&_r);
>>   VecDestroy(&_y);
>>   VecDestroy(&_z);
>>   VecDestroy(&_y0);
>>   VecDestroy(&_z0);
>>
>>   // set map iterator, free memory from VecScatter
>>   map<string,VecScatter>::iterator it;
>>   for (it = _scatters.begin(); it != _scatters.end(); it++) {
>>     VecScatterDestroy(&(it->second));
>>   }
>> }
>>
>> // construct coordinate transform, setting vectors q, r, y, z 
>> PetscErrorCode Domain::setFields() {
>>   PetscErrorCode ierr = 0;
>>
>>   // generate vector _y with size _Ny*_Nz
>>   ierr = VecCreate(PETSC_COMM_WORLD,&_y); CHKERRQ(ierr);
>>   ierr = VecSetSizes(_y,PETSC_DECIDE,_Ny*_Nz); CHKERRQ(ierr);
>>   ierr = VecSetFromOptions(_y); CHKERRQ(ierr);
>>   ierr = PetscObjectSetName((PetscObject) _y, "y"); CHKERRQ(ierr);
>>
>>   // duplicate _y into _z, _q, _r
>>   ierr = VecDuplicate(_y,&_z); CHKERRQ(ierr);
>>   ierr = PetscObjectSetName((PetscObject) _z, "z"); CHKERRQ(ierr);
>>   ierr = VecDuplicate(_y,&_q); CHKERRQ(ierr);
>>   ierr = PetscObjectSetName((PetscObject) _q, "q"); CHKERRQ(ierr);
>>   ierr = VecDuplicate(_y,&_r); CHKERRQ(ierr);
>>   ierr = PetscObjectSetName((PetscObject) _r, "r"); CHKERRQ(ierr);
>>
>>   // construct coordinate transform
>>   PetscInt Ii,Istart,Iend,Jj = 0;
>>   PetscScalar *y,*z,*q,*r;
>>   ierr = VecGetOwnershipRange(_q,&Istart,&Iend);CHKERRQ(ierr);
>>
>>   // return pointers to local data arrays (the processor's portion of vector data)
>>   ierr = VecGetArray(_y,&y); CHKERRQ(ierr);
>>   ierr = VecGetArray(_z,&z); CHKERRQ(ierr);
>>   ierr = VecGetArray(_q,&q); CHKERRQ(ierr);
>>   ierr = VecGetArray(_r,&r); CHKERRQ(ierr);
>>
>>   // set vector entries for q, r (coordinate transform) and y, z (no transform)
>>   for (Ii=Istart; Ii<Iend; Ii++) {
>>     q[Jj] = _dq*(Ii/_Nz);
>>     r[Jj] = _dr*(Ii-_Nz*(Ii/_Nz));
>>
>>     // matrix-based, fully compatible, allows curvilinear coordinate transformation
>>     if (_sbpType.compare("mfc_coordTrans") ) {
>>       y[Jj] = (_dq*_Ly)*(Ii/_Nz);
>>       z[Jj] = (_dr*_Lz)*(Ii-_Nz*(Ii/_Nz));
>>     }
>>     else {
>>       // hardcoded transformation (not available for z)
>>       if (_bCoordTrans > 0) {
>> 	y[Jj] = _Ly * sinh(_bCoordTrans * q[Jj]) / sinh(_bCoordTrans);
>>       }
>>       // no transformation
>>       y[Jj] = q[Jj]*_Ly;
>>       z[Jj] = r[Jj]*_Lz;
>>     }
>>     Jj++;
>>   }
>>
>>   // restore arrays
>>   ierr = VecRestoreArray(_y,&y); CHKERRQ(ierr);
>>   ierr = VecRestoreArray(_z,&z); CHKERRQ(ierr);
>>   ierr = VecRestoreArray(_q,&q); CHKERRQ(ierr);
>>   ierr = VecRestoreArray(_r,&r); CHKERRQ(ierr);
>>
>>   return ierr;
>> }
>>
>>
>> // scatters values from one vector to another PetscErrorCode 
>> Domain::setScatters() {
>>   PetscErrorCode ierr = 0;
>>
>>   ierr = VecCreate(PETSC_COMM_WORLD,&_y0); CHKERRQ(ierr);
>>   ierr = VecSetSizes(_y0,PETSC_DECIDE,_Nz); CHKERRQ(ierr);
>>   ierr = VecSetFromOptions(_y0); CHKERRQ(ierr);
>>   ierr = VecSet(_y0,0.0); CHKERRQ(ierr);
>>
>>   ierr = VecCreate(PETSC_COMM_WORLD,&_z0); CHKERRQ(ierr);
>>   ierr = VecSetSizes(_z0,PETSC_DECIDE,_Ny); CHKERRQ(ierr);
>>   ierr = VecSetFromOptions(_z0); CHKERRQ(ierr);
>>   ierr = VecSet(_z0,0.0); CHKERRQ(ierr);
>>
>>   PetscInt *indices;
>>   IS is;
>>   ierr = PetscMalloc1(_Nz,&indices); CHKERRQ(ierr);
>>
>>   // we want to scatter from index 0 to _Nz - 1, i.e. take the first _Nz components of the vector to scatter from
>>   for (PetscInt Ii = 0; Ii<_Nz; Ii++) {
>>     indices[Ii] = Ii;
>>   }
>>
>>   ierr = ISCreateGeneral(PETSC_COMM_WORLD, _Nz, indices, PETSC_COPY_VALUES, &is); CHKERRQ(ierr);
>>   ierr = VecScatterCreate(_y, is, _y0, is, &_scatters["body2L"]); 
>> CHKERRQ(ierr);
>>
>>   // free memory
>>   ierr = PetscFree(indices); CHKERRQ(ierr);
>>   ierr = ISDestroy(&is); CHKERRQ(ierr);
>>
>>   //===============================================================================
>>   // set up scatter context to take values for y = Ly from body field and put them on a Vec of size Nz
>>   PetscInt *fi;
>>   IS isf;
>>   ierr = PetscMalloc1(_Nz,&fi); CHKERRQ(ierr);
>>
>>   for (PetscInt Ii = 0; Ii<_Nz; Ii++) {
>>     fi[Ii] = Ii + (_Ny*_Nz-_Nz);
>>   }
>>   ierr = ISCreateGeneral(PETSC_COMM_WORLD, _Nz, fi, PETSC_COPY_VALUES, 
>> &isf); CHKERRQ(ierr);
>>
>>   PetscInt *ti;
>>   IS ist;
>>   ierr = PetscMalloc1(_Nz,&ti); CHKERRQ(ierr);
>>   for (PetscInt Ii = 0; Ii<_Nz; Ii++) {
>>     ti[Ii] = Ii;
>>   }
>>   ierr = ISCreateGeneral(PETSC_COMM_WORLD, _Nz, ti, PETSC_COPY_VALUES, &ist); CHKERRQ(ierr);
>>   ierr = VecScatterCreate(_y, isf, _y0, ist, &_scatters["body2R"]); 
>> CHKERRQ(ierr);
>>
>>   // free memory
>>   ierr = PetscFree(fi); CHKERRQ(ierr);
>>   ierr = PetscFree(ti); CHKERRQ(ierr);
>>   ierr = ISDestroy(&isf); CHKERRQ(ierr);
>>   ierr = ISDestroy(&ist); CHKERRQ(ierr);
>>
>>   
>>   //============================================================================== 
>>   IS isf2;
>>   ierr = ISCreateStride(PETSC_COMM_WORLD, _Ny, 0, _Nz, &isf2); 
>> CHKERRQ(ierr);
>>
>>   PetscInt *ti2;
>>   IS ist2;
>>   ierr = PetscMalloc1(_Ny,&ti2); CHKERRQ(ierr);
>>
>>   for (PetscInt Ii=0; Ii<_Ny; Ii++) {
>>     ti2[Ii] = Ii;
>>   }
>>   ierr = ISCreateGeneral(PETSC_COMM_WORLD, _Ny, ti2, PETSC_COPY_VALUES, &ist2); CHKERRQ(ierr);
>>   ierr = VecScatterCreate(_y, isf2, _z0, ist2, &_scatters["body2T"]); 
>> CHKERRQ(ierr);
>>
>>   // free memory
>>   ierr = PetscFree(ti2); CHKERRQ(ierr);
>>   ierr = ISDestroy(&isf2); CHKERRQ(ierr);
>>   ierr = ISDestroy(&ist2); CHKERRQ(ierr);
>>
>>
>>   //==============================================================================
>>   IS isf3;
>>   ierr = ISCreateStride(PETSC_COMM_WORLD, _Ny, _Nz - 1, _Nz, &isf3); 
>> CHKERRQ(ierr);
>>
>>   PetscInt *ti3;
>>   IS ist3;
>>   ierr = PetscMalloc1(_Ny,&ti3); CHKERRQ(ierr);
>>   for (PetscInt Ii = 0; Ii<_Ny; Ii++) {
>>     ti3[Ii] = Ii;
>>   }
>>   ierr = ISCreateGeneral(PETSC_COMM_WORLD, _Ny, ti3, PETSC_COPY_VALUES, &ist3); CHKERRQ(ierr);
>>   ierr = VecScatterCreate(_y, isf3, _z0, ist3, &_scatters["body2B"]); 
>> CHKERRQ(ierr);
>>
>>   // free memory
>>   ierr = PetscFree(ti3); CHKERRQ(ierr);
>>   ierr = ISDestroy(&isf3); CHKERRQ(ierr);
>>   ierr = ISDestroy(&ist3); CHKERRQ(ierr);
>>
>>   return ierr;
>> }
>> #include "domain.hpp"
>>
>> using namespace std;
>>
>> // creates a domain object
>> int main(int argc, char **argv) {
>>
>>   PetscErrorCode ierr = 0;
>>   PetscInitialize(&argc, &argv, NULL, NULL);
>>
>>   Domain d;
>>
>>   ierr = d.setFields(); CHKERRQ(ierr);
>>   ierr = d.setScatters(); CHKERRQ(ierr);
>>   
>>   PetscFinalize();
>>   return ierr;
>> }


More information about the petsc-users mailing list