[petsc-users] error: invalid types ‘PetscScalar* {aka double*}[PetscScalar {aka double}]’ for array subscript
Matthew Knepley
knepley at gmail.com
Fri Aug 21 11:19:34 CDT 2020
On Fri, Aug 21, 2020 at 11:49 AM Smit Thijs <thijs.smit at hest.ethz.ch> wrote:
> Hi All,
>
>
>
> I am having the following error when I try to do a mapping with vectors
> and I can’t figure out how to solve this or what is going wrong:
>
> error: invalid types ‘PetscScalar* {aka double*}[PetscScalar {aka
> double}]’ for array subscript
>
> xpMMA[i] = xp[indicesMap[i]];
>
>
>
> Herewith two code snippets:
>
> // total number of elements on core
>
> PetscInt nel;
>
> VecGetLocalSize(xPhys, &nel);
>
>
>
> // create xPassive vector
>
> ierr = VecDuplicate(xPhys, &xPassive);
>
> CHKERRQ(ierr);
>
>
>
> // create mapping vector
>
> ierr = VecDuplicate(xPhys, &indicator);
>
> CHKERRQ(ierr);
>
>
>
> // index set for xPassive and indicator
>
> PetscScalar *xpPassive, *xpIndicator;
>
> ierr = VecGetArray(xPassive, &xpPassive);
>
> CHKERRQ(ierr);
>
> ierr = VecGetArray(indicator, &xpIndicator);
>
> CHKERRQ(ierr);
>
>
>
> // counters for total and active elements on this processor
>
> PetscInt tcount = 0; // total number of elements
>
> PetscInt acount = 0; // number of active elements
>
> PetscInt scount = 0; // number of solid elements
>
> PetscInt rcount = 0; // number of rigid element
>
>
>
> // loop over all elements and update xPassive from wrapper data
>
> // count number of active elements, acount
>
> // set indicator vector
>
> for (PetscInt el = 0; el < nel; el++) {
>
> if (data.xPassive_w.size() > 1) {
>
> xpPassive[el] = data.xPassive_w[el];
>
> tcount++;
>
> if (xpPassive[el] < 0) {
>
> xpIndicator[acount] = el;
>
> acount++;
>
> }
>
> } else {
>
> xpPassive[el] = -1.0; // default, if no xPassive_w than all
> elements are active = -1.0
>
> }
>
> }
>
>
>
> // printing
>
> //PetscPrintf(PETSC_COMM_WORLD, "tcount: %i\n", tcount);
>
> //PetscPrintf(PETSC_COMM_WORLD, "acount: %i\n", acount);
>
>
>
> // Allreduce, get number of active elements over all processes
>
> // tmp number of var on proces
>
> // acount total number of var sumed
>
> PetscInt tmp = acount;
>
> acount = 0.0;
>
> MPI_Allreduce(&tmp, &(acount), 1, MPIU_INT, MPI_SUM, PETSC_COMM_WORLD);
>
>
>
> //// create xMMA vector
>
> VecCreateMPI(PETSC_COMM_WORLD, tmp, acount, &xMMA);
>
>
>
> // Pointers to the vectors
>
> PetscScalar *xp, *xpMMA, *indicesMap;
>
Here you declare indicesMap as PetscScalar[]. You cannot index an array
with this. I see that you
want to store these indices in a Vec. You should use an IS instead.
Thanks,
Matt
> //PetscInt indicesMap;
>
> ierr = VecGetArray(MMAVector, &xpMMA);
>
> CHKERRQ(ierr);
>
> ierr = VecGetArray(elementVector, &xp);
>
> CHKERRQ(ierr);
>
> // Index set
>
> PetscInt nLocalVar;
>
> VecGetLocalSize(xMMA, &nLocalVar);
>
>
>
> // print number of var on pocessor
>
> PetscPrintf(PETSC_COMM_WORLD, "Local var: %i\n", nLocalVar);
>
>
>
> ierr = VecGetArray(indicator, &indicesMap);
>
> CHKERRQ(ierr);
>
>
>
> // Run through the indices
>
> for (PetscInt i = 0; i < nLocalVar; i++) {
>
> if (updateDirection > 0) {
>
> //PetscPrintf(PETSC_COMM_WORLD, "i: %i, xp[%i] = %f\n", i,
> indicesMap[i], xp[indicesMap[i]]);
>
> xpMMA[i] = xp[indicesMap[i]];
>
> } else if (updateDirection < 0) {
>
> xp[indicesMap[i]] = xpMMA[i];
>
> //PetscPrintf(PETSC_COMM_WORLD, "i: %i, xp[%i] = %f\n", i,
> indicesMap[i], xp[indicesMap[i]]);
>
> }
>
> }
>
> // Restore
>
> ierr = VecRestoreArray(elementVector, &xp);
>
> CHKERRQ(ierr);
>
> ierr = VecRestoreArray(MMAVector, &xpMMA);
>
> CHKERRQ(ierr);
>
> ierr = VecRestoreArray(indicator, &indicesMap);
>
> CHKERRQ(ierr);
>
> PetscPrintf(PETSC_COMM_WORLD, "FINISHED UpdateVariables \n");
>
>
>
> The error message says that the type with which I try to index is wrong, I
> think. But VecGetArray only excepts scalars. Furthermore, the el variable
> is an int, but is seams like to turn out to be a scalar. Does anybody see
> how to proceed with this?
>
>
>
> Best regards,
>
>
>
> Thijs Smit
>
--
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/~knepley/>
-------------- next part --------------
An HTML attachment was scrubbed...
URL: <http://lists.mcs.anl.gov/pipermail/petsc-users/attachments/20200821/f8d4a19b/attachment.html>
More information about the petsc-users
mailing list