Hi, PETSc developer:

I need to create a PetscSection with a struct. The struct is defined as follow,

   typedef struct
   {
      PetscReal x;
      PetscInt id;
   } testStruct;

   When I run the program, I got a wrong output as follow,

   Vec Object: 1 MPI processes
  type: seq
2.
4.94066e-324
2.
4.94066e-324
2.
4.94066e-324
2.
4.94066e-324
2.
4.94066e-324
2.
4.94066e-324
2.
4.94066e-324
2.
4.94066e-324

But when I defined the struct as

   typedef struct
   {
      PetscReal x;
      PetscReal id;
   } testStruct;

The output is ok. It seems that there is some wrong with the memories when I define the "id" as a PetscInt type.

I can not find out the reasons, and any one can help me with it? The source file "test.c" is attached.


Thanks,

leejearl

static char help[] = "An example \n\n";

#include <petscdmplex.h>

typedef struct
{
    PetscReal x;
    PetscInt id;
} testStruct;

int main(int argc, char **argv)
{
  DM             dm, dmDist = NULL;
  Vec            u;
  PetscSection   section;
  PetscViewer    viewer;
  PetscInt       dim = 2, numFields, numBC, i;
  PetscInt       numComp[3];
  PetscInt       numDof[12];
  PetscInt       bcField[1];
  IS             bcPointIS[1];
  PetscBool      interpolate = PETSC_TRUE;
  PetscErrorCode ierr;

  ierr = PetscInitialize(&argc, &argv, NULL,help);if (ierr) return ierr;
  /* Create a mesh */
  ierr = DMPlexCreateBoxMesh(PETSC_COMM_WORLD, dim, 2, interpolate, &dm);CHKERRQ(ierr);

  {
      PetscInt c, cStart, cEnd, cEndInterior;
      Vec testVec;
      PetscScalar *testScalar;
      DM dmTest;
      PetscSection   coordSection, sectionTest;
      Vec            coordinates;
      DMGetDimension(dm, &dim);
      DMGetCoordinateSection(dm, &coordSection);
      DMGetCoordinatesLocal(dm, &coordinates);

      DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);
      DMPlexGetHybridBounds(dm, &cEndInterior, NULL, NULL, NULL);

      DMClone(dm, &dmTest);
      DMSetCoordinateSection(dmTest, PETSC_DETERMINE, coordSection);
      DMSetCoordinatesLocal(dmTest, coordinates);
      PetscSectionCreate(PetscObjectComm((PetscObject) dm), &sectionTest);
      PetscSectionSetChart(sectionTest, cStart, cEnd);
      for (c = cStart; c < cEnd; ++c) 
      {
          PetscSectionSetDof(sectionTest, c, 2);
      }
      PetscSectionSetUp(sectionTest);
      DMSetDefaultSection(dmTest, sectionTest);
      PetscSectionDestroy(&sectionTest);

      DMCreateLocalVector(dmTest, &testVec);
      if (cEndInterior < 0) {
          cEndInterior = cEnd;
      }
      VecGetArray(testVec, &testScalar);
      for (c = cStart; c < cEndInterior; ++c) {
          testStruct *testData;

          DMPlexPointLocalRef(dmTest, c, testScalar, &testData);
          testData->x = 2.0;
          testData->id = 1;
      }
      VecRestoreArray(testVec, &testScalar);

      VecView(testVec,PETSC_VIEWER_STDOUT_WORLD);
  }

 
  ierr = DMDestroy(&dm);CHKERRQ(ierr);
  ierr = PetscFinalize();
  return ierr;
}

Reply via email to