[petsc-users] Issue with Exporting and Reading Complex Vectors and Magnitude Vectors in PETSc and MATLAB

Barry Smith bsmith at petsc.dev
Tue Jul 2 09:24:20 CDT 2024


arecomplex = false;

tnargin = nargin;
for l=1:nargin-2
  if ischar(varargin{l}) && strcmpi(varargin{l},'indices')
    tnargin = min(l,tnargin-1);
    indices = varargin{l+1};
  end
  if ischar(varargin{l}) && strcmpi(varargin{l},'precision')
    tnargin = min(l,tnargin-1);
    precision = varargin{l+1};
  end
  if ischar(varargin{l}) && strcmpi(varargin{l},'cell')
    tnargin = min(l,tnargin-1);
    arecell = varargin{l+1};
  end
  if ischar(varargin{l}) && strcmpi(varargin{l},'complex')      <======== finds any argument that is 'complex'
    tnargin = min(l,tnargin-1);
    arecomplex = varargin{l+1};                                            <======== sets arecomplex to the next argument in the list
  end
end

If you are still having trouble you can use the Matlab debugger to step through this code to see why this check is not triggered.



> On Jul 2, 2024, at 4:59 AM, maitri ksh <maitri.ksh at gmail.com> wrote:
> 
> This Message Is From an External Sender
> This message came from outside your organization.
> A variable 'arecomplex' which is by default set to false inside PetscBinaryRead() overrides 'complex', true in PetscBinaryRead(file, "complex", true) thus giving real valued output unless one manually changes arecomplex=true  inside the function.
> 
> On Mon, Jul 1, 2024 at 5:40 PM Pierre Jolivet <pierre at joliv.et <mailto:pierre at joliv.et>> wrote:
>> 
>> 
>>> On 1 Jul 2024, at 4:37 PM, maitri ksh <maitri.ksh at gmail.com <mailto:maitri.ksh at gmail.com>> wrote:
>>> 
>>> This Message Is From an External Sender
>>> This message came from outside your organization.
>>> I need to export complex vectors data from PETSc and read them in MATLAB. However, I am encountering some issues with the data format and interpretation in MATLAB.  
>>> 
>>> code-snippet of the vector data export section:
>>> // Assemble the vectors before exporting
>>> ierr = VecAssemblyBegin(f); CHKERRQ(ierr);
>>> ierr = VecAssemblyEnd(f); CHKERRQ(ierr);
>>> 
>>> PetscViewer viewerf;
>>> // Save the complex vectors to binary files
>>> ierr = PetscViewerBinaryOpen(PETSC_COMM_WORLD, "output_f.dat", FILE_MODE_WRITE, &viewerf); CHKERRQ(ierr);
>>> ierr = VecView(f, viewerf); CHKERRQ(ierr);
>>> ierr = PetscViewerDestroy(&viewerf); CHKERRQ(ierr);
>>> 
>>> // Create vectors to store the magnitudes
>>> Vec f_magnitude;
>>> ierr = VecDuplicate(f, &f_magnitude); CHKERRQ(ierr);
>>> 
>>> // Get local portion of the vectors
>>> const PetscScalar *f_array;
>>> PetscScalar *f_magnitude_array;
>>> PetscInt n_local;
>>> 
>>> ierr = VecGetLocalSize(f, &n_local); CHKERRQ(ierr);
>>> ierr = VecGetArrayRead(f, &f_array); CHKERRQ(ierr);
>>> ierr = VecGetArray(f_magnitude, &f_magnitude_array); CHKERRQ(ierr);
>>> 
>>> // Compute the magnitude for each element
>>> for (int i = 0; i < n_local; i++) {
>>>     f_magnitude_array[i] = PetscAbsScalar(f_array[i]);
>>> }
>>> 
>>> // Restore arrays
>>> ierr = VecRestoreArrayRead(f, &f_array); CHKERRQ(ierr);
>>> ierr = VecRestoreArray(f_magnitude, &f_magnitude_array); CHKERRQ(ierr);
>>> 
>>> // Assemble the magnitude vectors
>>> ierr = VecAssemblyBegin(f_magnitude); CHKERRQ(ierr);
>>> ierr = VecAssemblyEnd(f_magnitude); CHKERRQ(ierr);
>>> 
>>> // Save the magnitude vectors to binary files
>>> PetscViewer viewerfmag;
>>> ierr = PetscViewerBinaryOpen(PETSC_COMM_WORLD, "output_f_mag.dat", FILE_MODE_WRITE, &viewerfmag); CHKERRQ(ierr);
>>> ierr = VecView(f_magnitude, viewerfmag); CHKERRQ(ierr);
>>> ierr = PetscViewerDestroy(&viewerfmag); CHKERRQ(ierr);
>>> 
>>> In MATLAB, I am using petscBinaryRead to read the data. The complex vectors are read, but only the real part is available. The magnitude vectors are however read as alternating zero and non-zero elements. What went wrong?
>>> How can I export the data correctly to MATLAB-accessible format? (I have not configured PETSc with Matlab as I was encountering library conflict issues)
>> 
>> 
>> How are you reading complex-valued binary files in MATLAB?
>> It is not the same as real-valued files, in particular, you must add the parameters 'complex', true to the PetscBinaryRead() call, type help PetscBinaryRead in the MATLAB console.
>> 
>> Thanks,
>> Pierre
>> 
> 

-------------- next part --------------
An HTML attachment was scrubbed...
URL: <http://lists.mcs.anl.gov/pipermail/petsc-users/attachments/20240702/3cf0d4c2/attachment-0001.html>


More information about the petsc-users mailing list