[petsc-users] valgrind errors.

(Rebecca) Xuefei YUAN xy2102 at columbia.edu
Tue Jan 19 12:49:36 CST 2010


Dear Barry,

I switched from using temp arrays via PetscMalloc() and PetscFree() to  
creating a DA and using DAGetLocalVector()+DAVecGetArray() to get  
fieldother as a temp array.

However, I do not think DAGetLocalVector()+DAVecGetArray() is a must  
in my case, I tried DAGetArray() to take over of the above two.  
However, I have the error message saying


twgcqt2unff.c:2150: warning: passing argument 3 of ?DAGetArray? from  
incompatible pointer type
twgcqt2unff.c:3213: warning: passing argument 3 of ?DARestoreArray?  
from incompatible pointer type
.

The code looks like:



typedef struct {
	PetscReal	x1,x2,x3,x4;
} FieldOther;

#undef __FUNCT__
#define __FUNCT__ "FormFunction"
PetscErrorCode FormFunction(SNES snes,Vec X,Vec F,void*dummg)
{
// ...
	DA				dafieldother;
	FieldOther		**fieldother;
	Vec				FIELDOTHER;
// ...
	PetscFunctionBegin;
// ...
	ierr = DAGetLocalInfo((DA)dmmg->dm,&info);CHKERRQ(ierr);
// ...
	ierr = DACreate2d(PETSC_COMM_WORLD,DA_NONPERIODIC,DA_STENCIL_BOX,  
info.mx, info.my, PETSC_DECIDE, PETSC_DECIDE, 4, 2, 0, 0,  
&dafieldother);CHKERRQ(ierr);
	ierr = DASetFieldName(dafieldother, 0, "x1");CHKERRQ(ierr);
	ierr = DASetFieldName(dafieldother, 1, "x2");CHKERRQ(ierr);
	ierr = DASetFieldName(dafieldother, 2, "x3");CHKERRQ(ierr);
	ierr = DASetFieldName(dafieldother, 3, "x4");CHKERRQ(ierr);

//	ierr = DAGetLocalVector(dafieldother,&FIELDOTHER);CHKERRQ(ierr);
//  ierr = DAVecGetArray(dafieldother,FIELDOTHER,&fieldother);CHKERRQ(ierr);

	ierr = DAGetArray(dafieldother,PETSC_TRUTH,&fieldother);CHKERRQ(ierr);

	for (j = jFirst; j <= jLast; j++){
         for (i = iFirst; i <= iLast; i++){
                 fieldother[j][i].x1 = 0;
                 fieldother[j][i].x2 = 0;
                 fieldother[j][i].x3 = 0;
                 fieldother[j][i].x4 = 0;
		}
	}
// ...
     for (j = jFirst; j <= jLast; j++){
         for (i = iFirst; i <= iLast; i++){
			ffield[j][i].x1 = fieldother[j][i].x1;
		}
	}
// ...
//  ierr =  
DAVecRestoreArray(dafieldother,FIELDOTHER,&fieldother);CHKERRQ(ierr);
//	ierr = DARestoreLocalVector(dafieldother,&FIELDOTHER);CHKERRQ(ierr);

	ierr = DARestoreArray(dafieldother,PETSC_TRUTH,&fieldother);CHKERRQ(ierr);

	ierr = DADestroy(dafieldother);CHKERRQ(ierr);
	PetscFunctionReturn(0);

}

What could be wrong with it?

Thanks very much!

Rebecca

Quoting Barry Smith <bsmith at mcs.anl.gov>:

>
>> err =   
>> PetscMalloc(sizeof(PetscReal)*(info.mx)*(info.my)*sizeof(FieldOther),   
>> &fieldother);CHKERRQ(ierr);
>
>     This line is wrong, but is not the cause of the problem.
> sizeof(FieldOther) takes into account the size of the real values
> inside so you do not need the sizeof(PetscReal)*
>
>   What are iFirst and iLast? If they are not 0 then you indexing into
> fieldother[] wrong since that array starts at 0.
>
>    Barry
>
> On Jan 19, 2010, at 11:24 AM, (Rebecca) Xuefei YUAN wrote:
>
>> Dear Jed,
>>
>> I found the source of these uninitialised values.
>>
>> In my FormFunction(), there is a user-defined structure called   
>> FieldOther and I will use this as a temp values to build up my   
>> residual function. The main loop is like:
>>
>>
>> typedef struct {
>> 	PetscReal	x1,x2,x3,x4;
>> } FieldOther;
>>
>> #undef __FUNCT__
>> #define __FUNCT__ "FormFunction"
>> PetscErrorCode FormFunction(SNES snes,Vec X,Vec F,void*dummg)
>> {
>> // ...
>> 	FieldOther		*fieldother;
>> // ...
>> 	PetscFunctionBegin;
>> // ...
>> 	ierr =   
>> PetscMalloc(sizeof(PetscReal)*(info.mx)*(info.my)*sizeof(FieldOther),   
>> &fieldother);CHKERRQ(ierr);
>> // ...
>>   for (j = jFirst; j <= jLast; j++){
>>       for (i = iFirst; i <= iLast; i++){
>>               fieldother[(i)+(j)*(info.mx)].x1 = 0;
>>               fieldother[(i)+(j)*(info.mx)].x2 = 0;
>>               fieldother[(i)+(j)*(info.mx)].x3 = 0;
>>               fieldother[(i)+(j)*(info.mx)].x4 = 0;
>> 		}
>> 	}
>> // ...
>>   for (j = jFirst; j <= jLast; j++){
>>       for (i = iFirst; i <= iLast; i++){
>> // errors coming out from valgrind if fieldother[] is used in   
>> residual function evaluation.
>> //			ffield[j][i].x1 = fieldother[(i)+(j)*(info.mx)].x1;
>> // everything looks fine from valgrind if fieldother[] is not used   
>> in residual function evaluation.
>> 			ffield[j][i].x1 = field[j][i].x1;
>> 		}
>> 	}
>> // ...
>> 	ierr = PetscFree(fieldother);CHKERRQ(ierr);
>> // ...
>> 	PetscFunctionReturn(0);
>>
>> }
>>
>> If I will use those fieldother[] in my residual evaluation   
>> functions, errors from valgrind come out like:
>>
>> ==24486==
>> ==24486== Conditional jump or move depends on uninitialised value(s)
>> ==24486==    at 0x4024D9B: memcpy (mc_replace_strmem.c:482)
>> ==24486==    by 0x41F6387: __printf_fp (in /lib/tls/i686/cmov/libc-2.7.so)
>> ==24486==    by 0x41F1563: vfprintf (in /lib/tls/i686/cmov/libc-2.7.so)
>> ==24486==    by 0x41F5581: ??? (in /lib/tls/i686/cmov/libc-2.7.so)
>> ==24486==    by 0x41F1229: vfprintf (in /lib/tls/i686/cmov/libc-2.7.so)
>> ==24486==    by 0x8643996: PetscVFPrintfDefault (mprint.c:193)
>> ==24486==    by 0x85F3F00: PetscViewerASCIIMonitorPrintf (filev.c:988)
>> ==24486==    by 0x80C814A: SNESMonitorDefault (snesut.c:150)
>> ==24486==    by 0x80E0450: SNESSolve_LS (ls.c:173)
>> ==24486==    by 0x80C199F: SNESSolve (snes.c:2221)
>> ==24486==    by 0x80DA4A0: DMMGSolveSNES (damgsnes.c:510)
>> ==24486==    by 0x80D3CDF: DMMGSolve (damg.c:372)
>> ==24486==  Uninitialised value was created by a heap allocation
>> ==24486==    at 0x4022BF3: malloc (vg_replace_malloc.c:195)
>> ==24486==    by 0x864BAFC: PetscMallocAlign (mal.c:40)
>> ==24486==    by 0x864CBEB: PetscTrMallocDefault (mtr.c:194)
>> ==24486==    by 0x8055FB8: FormFunction (twgcqt2unff.c:941)
>> ==24486==    by 0x80B86A4: SNESComputeFunction (snes.c:1016)
>> ==24486==    by 0x80DFFC2: SNESSolve_LS (ls.c:159)
>> ==24486==    by 0x80C199F: SNESSolve (snes.c:2221)
>> ==24486==    by 0x80DA4A0: DMMGSolveSNES (damgsnes.c:510)
>> ==24486==    by 0x80D3CDF: DMMGSolve (damg.c:372)
>> ==24486==    by 0x804F539: Solve (twgcqt2unff.c:437)
>> ==24486==    by 0x804DC92: main (twgcqt2unff.c:283)
>>
>> FormFunction (twgcqt2unff.c:941) is located at the line of
>>
>> 	ierr =   
>> PetscMalloc(sizeof(PetscReal)*(info.mx)*(info.my)*sizeof(FieldOther),   
>> &fieldother);CHKERRQ(ierr);
>>
>> I think I did give values of those fieldother and free them after   
>> usage, but why this uninitialised message comes out all the time?
>>
>> Thanks very much!
>>
>> Rebecca
>>
>>
>> Quoting Jed Brown <jed at 59A2.org>:
>>
>>> On Mon, 18 Jan 2010 16:32:22 -0500, "(Rebecca) Xuefei YUAN"    
>>> <xy2102 at columbia.edu> wrote:
>>>> Dear Jed,
>>>>
>>>> I got some origins, but I still do not understand what is going on here?
>>>>
>>>> Here is a piece of uninitialised value from
>>>
>>> I think you have an indexing error while reading from this local vector
>>> in your residual evaluation.  The contents of the array are explicitly
>>> initialized to zero in VecCreate_Seq, but when memalign is not available
>>> (your trace indicates this is so) PetscMallocAlign allocate a little
>>> extra space and aligns the result, leaving a few bytes uninitialized
>>> (memory further away may also be uninitialized, but this error isn't
>>> saying complaining about those places).
>>>
>>> Jed
>>>
>>>
>>
>>
>>
>> -- 
>> (Rebecca) Xuefei YUAN
>> Department of Applied Physics and Applied Mathematics
>> Columbia University
>> Tel:917-399-8032
>> www.columbia.edu/~xy2102
>>



-- 
(Rebecca) Xuefei YUAN
Department of Applied Physics and Applied Mathematics
Columbia University
Tel:917-399-8032
www.columbia.edu/~xy2102



More information about the petsc-users mailing list