<html>
  <head>
    <meta http-equiv="Content-Type" content="text/html; charset=utf-8">
  </head>
  <body text="#000000" bgcolor="#FFFFFF">
    On 18-02-20 09:52 AM, Matthew Knepley wrote:<br>
    <blockquote type="cite"
cite="mid:CAMYG4GknBNzGcKs-SuUwp_CEVPXB7n65DZbOrLWFjjA9iMTx9g@mail.gmail.com">
      <div dir="ltr">
        <div class="gmail_extra">
          <div class="gmail_quote">On Tue, Feb 20, 2018 at 12:30 PM,
            Danyang Su <span dir="ltr"><<a
                href="mailto:danyang.su@gmail.com" target="_blank"
                moz-do-not-send="true">danyang.su@gmail.com</a>></span>
            wrote:<br>
            <blockquote class="gmail_quote" style="margin:0px 0px 0px
              0.8ex;border-left:1px solid
              rgb(204,204,204);padding-left:1ex">
              <div bgcolor="#FFFFFF">
                <p>Hi All,</p>
                <p>I tried to compile the DMPlexCreateSection code but
                  got error information as shown below.<br>
                </p>
                <p>Error: Symbol 'petsc_null_is' at (1) has no IMPLICIT
                  type</p>
                <p>I tried to use PETSC_NULL_OBJECT instead of
                  PETSC_NULL_IS, then the code can be compiled but run
                  into Segmentation Violation error in
                  DMPlexCreateSection.<br>
                </p>
              </div>
            </blockquote>
            <div>From the webpage</div>
            <div><br>
            </div>
            <div>  <a
href="http://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/DMPLEX/DMPlexCreateSection.html"
                moz-do-not-send="true">http://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/DMPLEX/DMPlexCreateSection.html</a> <br>
            </div>
            <div><br>
            </div>
            <div>The F90 version is DMPlexCreateSectionF90. Doing this
              with F77 arrays would have been too painful.</div>
          </div>
        </div>
      </div>
    </blockquote>
    Hi Matt,<br>
    <br>
    Sorry, I still cannot compile the code if use DMPlexCreateSectionF90
    instead of DMPlexCreateSection. Would you please tell me in more
    details? <br>
    <br>
    undefined reference to `dmplexcreatesectionf90_'<br>
    <br>
    then I #include <petsc/finclude/petscdmplex.h90>, but this
    throws more error during compilation.<br>
    <br>
    <br>
        Included at
/home/dsu/Soft/PETSc/petsc-3.7.5/include/petsc/finclude/petscdmplex.h90:6:<br>
        Included at ../../solver/solver_ddmethod.F90:62:<br>
    <br>
              PETSCSECTION_HIDE section<br>
              1<br>
    Error: Unclassifiable statement at (1)<br>
/home/dsu/Soft/PETSc/petsc-3.7.5/include/petsc/finclude/ftn-custom/petscdmplex.h90:167.10:<br>
        Included at
/home/dsu/Soft/PETSc/petsc-3.7.5/include/petsc/finclude/petscdmplex.h90:6:<br>
        Included at ../../solver/solver_ddmethod.F90:62:<br>
    <br>
              PETSCSECTION_HIDE section<br>
              1<br>
    Error: Unclassifiable statement at (1)<br>
/home/dsu/Soft/PETSc/petsc-3.7.5/include/petsc/finclude/ftn-custom/petscdmplex.h90:179.10:<br>
        Included at
/home/dsu/Soft/PETSc/petsc-3.7.5/include/petsc/finclude/petscdmplex.h90:6:<br>
        Included at ../../solver/solver_ddmethod.F90:62:<br>
    <br>
    <blockquote type="cite"
cite="mid:CAMYG4GknBNzGcKs-SuUwp_CEVPXB7n65DZbOrLWFjjA9iMTx9g@mail.gmail.com">
      <div dir="ltr">
        <div class="gmail_extra">
          <div class="gmail_quote">
            <div><br>
            </div>
            <div>  Thanks,</div>
            <div><br>
            </div>
            <div>     Matt</div>
            <blockquote class="gmail_quote" style="margin:0px 0px 0px
              0.8ex;border-left:1px solid
              rgb(204,204,204);padding-left:1ex">
              <div bgcolor="#FFFFFF">
                <p> </p>
                <p>dmda_flow%da is distributed dm object that works
                  fine. </p>
                <p>The fortran example I follow is
                  <a
                    class="gmail-m_8962138310016109093moz-txt-link-freetext"
href="http://www.mcs.anl.gov/petsc/petsc-dev/src/dm/impls/plex/examples/tutorials/ex1f90.F90"
                    target="_blank" moz-do-not-send="true">http://www.mcs.anl.gov/petsc/<wbr>petsc-dev/src/dm/impls/plex/<wbr>examples/tutorials/ex1f90.F90</a>.
                  <br>
                </p>
                <p>What parameters should I use if passing null to
                  bcField, bcComps, bcPoints and perm.<br>
                </p>
                <pre style="color:rgb(0,0,0);font-style:normal;font-variant-ligatures:normal;font-variant-caps:normal;font-weight:400;letter-spacing:normal;text-align:start;text-indent:0px;text-transform:none;word-spacing:0px;text-decoration-style:initial;text-decoration-color:initial"><a href="https://www.mcs.anl.gov/petsc/petsc-master/docs/manualpages/Vec/PetscErrorCode.html#PetscErrorCode" target="_blank" moz-do-not-send="true">PetscErrorCode</a> <a href="https://www.mcs.anl.gov/petsc/petsc-master/docs/manualpages/DMPLEX/DMPlexCreateSection.html#DMPlexCreateSection" target="_blank" moz-do-not-send="true">DMPlexCreateSection</a>(<a href="https://www.mcs.anl.gov/petsc/petsc-master/docs/manualpages/DM/DM.html#DM" target="_blank" moz-do-not-send="true">DM</a> dm, <a href="https://www.mcs.anl.gov/petsc/petsc-master/docs/manualpages/Vec/PetscInt.html#PetscInt" target="_blank" moz-do-not-send="true">PetscInt</a> dim, <a href="https://www.mcs.anl.gov/petsc/petsc-master/docs/manualpages/Vec/PetscInt.html#PetscInt" target="_blank" moz-do-not-send="true">PetscInt</a> numFields,const <a href="https://www.mcs.anl.gov/petsc/petsc-master/docs/manualpages/Vec/PetscInt.html#PetscInt" target="_blank" moz-do-not-send="true">PetscInt</a> numComp[],const <a href="https://www.mcs.anl.gov/petsc/petsc-master/docs/manualpages/Vec/PetscInt.html#PetscInt" target="_blank" moz-do-not-send="true">PetscInt</a> numDof[], <a href="https://www.mcs.anl.gov/petsc/petsc-master/docs/manualpages/Vec/PetscInt.html#PetscInt" target="_blank" moz-do-not-send="true">PetscInt</a> numBC,const <a href="https://www.mcs.anl.gov/petsc/petsc-master/docs/manualpages/Vec/PetscInt.html#PetscInt" target="_blank" moz-do-not-send="true">PetscInt</a> bcField[], 
const <a href="https://www.mcs.anl.gov/petsc/petsc-master/docs/manualpages/IS/IS.html#IS" target="_blank" moz-do-not-send="true">IS</a> bcComps[], const <a href="https://www.mcs.anl.gov/petsc/petsc-master/docs/manualpages/IS/IS.html#IS" target="_blank" moz-do-not-send="true">IS</a> bcPoints[], <a href="https://www.mcs.anl.gov/petsc/petsc-master/docs/manualpages/IS/IS.html#IS" target="_blank" moz-do-not-send="true">IS</a> perm, <a href="https://www.mcs.anl.gov/petsc/petsc-master/docs/manualpages/IS/PetscSection.html#PetscSection" target="_blank" moz-do-not-send="true">PetscSection</a> *section)

</pre>
                <p>#include <petsc/finclude/petscis.h><br>
                  #include <petsc/finclude/petscis.h90><br>
                  #include <petsc/finclude/petscdmplex.h></p>
                <p>...<br>
                </p>
                <p>#ifdef USG<br>
                          numFields = 1<br>
                          numComp(1) = 1<br>
                          pNumComp => numComp<br>
                  <br>
                          do i = 1, numFields*(dmda_flow%dim+1)<br>
                            numDof(i) = 0<br>
                          end do<br>
                          numDof(0*(dmda_flow%dim+1)+1) = dmda_flow%dof<br>
                          pNumDof => numDof<br>
                  <br>
                          numBC = 0<br>
                  <br>
                          call DMPlexCreateSection(dmda_flow%<wbr>da,dmda_flow%dim,         
                  &<br>
                                                <wbr>  
                  numFields,pNumComp,pNumDof,   <wbr>                       
                  &<br>
                                                <wbr>  
                  numBC,PETSC_NULL_INTEGER,     <wbr>                         
                  &<br>
                                                <wbr>  
                  PETSC_NULL_IS,PETSC_NULL_IS,  <wbr>                          
                  &             !Error here<br>
                                                <wbr>  
                  PETSC_NULL_IS,section,ierr)<br>
                          CHKERRQ(ierr)<br>
                  <br>
                          call PetscSectionSetFieldName(<wbr>section,0,'flow',ierr)<br>
                          CHKERRQ(ierr)<br>
                  <br>
                          call DMSetDefaultSection(dmda_flow%<wbr>da,section,ierr)<br>
                          CHKERRQ(ierr)<br>
                  <br>
                          call PetscSectionDestroy(section,<wbr>ierr)<br>
                          CHKERRQ(ierr)<br>
                  #endif</p>
                <p>Thanks,</p>
                <p>Danyang<br>
                </p>
              </div>
            </blockquote>
          </div>
          <br>
          <br clear="all">
          <div><br>
          </div>
          -- <br>
          <div class="gmail_signature">
            <div dir="ltr">
              <div>
                <div dir="ltr">
                  <div>What most experimenters take for granted before
                    they begin their experiments is infinitely more
                    interesting than any results to which their
                    experiments lead.<br>
                    -- Norbert Wiener</div>
                  <div><br>
                  </div>
                  <div><a href="http://www.caam.rice.edu/%7Emk51/"
                      target="_blank" moz-do-not-send="true">https://www.cse.buffalo.edu/~knepley/</a><br>
                  </div>
                </div>
              </div>
            </div>
          </div>
        </div>
      </div>
    </blockquote>
    <br>
  </body>
</html>