[petsc-users] Compiling program with program files dependencies (SLEPC )

Steve Ndengue steve.ndengue at gmail.com
Fri Apr 25 09:45:57 CDT 2014


Dear all,

Sorry, I had a typo: I was still editing the .f90 program file instead 
of the .F90 file.

However, now I am having troubles forming the matrice I would like to 
diagonalize.
Presently the size is not yet huge and it could be held in a previously 
computed matrix of dimension N*N.

I am using the following to build the matrix and the compiling from the 
example 1 and slightly modified:

***
/!//
//        CALL SlepcInitialize(PETSC_NULL_CHARACTER,ierr)//
//        CALL MPI_Comm_rank(PETSC_COMM_WORLD,rank,ierr)//
//        nn = N//
//!        CALL PetscOptionsGetInt(PETSC_NULL_CHARACTER,'-n',n,flg,ierr)//
//!//
////
//!//
//! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 
- - //
//!     Compute the operator matrix that defines the eigensystem, Ax=kx//
//! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 
- - //
//
//      call MatCreate(PETSC_COMM_WORLD,A,ierr)//
//      call MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,nn,nn,ierr)//
//      call MatSetFromOptions(A,ierr)//
//      call MatSetUp(A,ierr)//
//
//      call MatGetOwnershipRange(A,Istart,Iend,ierr)//
//      do i=Istart,Iend-1//
//      do j=0,nn-1//
//      call MatSetValue(A,i,-1,j,-1,H_mat(i,j),INSERT_VALUES,ierr)//
//      enddo//
//      enddo//
//      call MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY,ierr)//
//      call MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY,ierr)! ** Create 
eigensolver context//
//      CALL EPSCreate(PETSC_COMM_WORLD,solver,ierr)//
//!//     ** Set operators. In this case, it is a standard eigenvalue 
problem//
//      CALL EPSSetOperators(solver,A,PETSC_NULL_OBJECT,ierr)//
//print*,'here10,ierr',ierr//
//      CALL EPSSetProblemType(solver,EPS_HEP,ierr)//
//!     ** Set working dimensions//
//      CALL EPSSETDimensions(solver,nstates,10*nstates,2*nstates)//
//!     ** Set calculation of lowest eigenvalues//
//      CALL EPSSetWhichEigenpairs(solver,EPS_SMALLEST_REAL)//
//
//!     ** Set solver parameters at runtime//
//      CALL EPSSetFromOptions(solver,ierr)//
//! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 
- - //
//!     Solve the eigensystem//
//! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 
- - //
//      CALL EPSSolve(solver,ierr)//
//print*,'here6'//
//!     ** Optional: Get some information from the solver and display it//
//      CALL EPSGetType(solver,tname,ierr)//
//      if (rank .eq. 0) then//
//        write(*,120) tname//
//      endif//
// 120  format (' Solution method: ',A)//
//      CALL EPSGetDimensions(solver,nstates,PETSC_NULL_INTEGER, &//
//     &                      PETSC_NULL_INTEGER,ierr)//
//      if (rank .eq. 0) then//
//        write(*,130) nstates//
//      endif//
// 130  format (' Number of requested eigenvalues:',I4)//
//! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 
- - //
//!     Display solution and clean up//
//! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 
- - //
//      CALL EPSPrintSolution(solver,PETSC_NULL_OBJECT,ierr)//
//      CALL EPSDestroy(solver,ierr)//
//***/

I get the following error message when making the code:

***
/*[0]PETSC ERROR: 
------------------------------------------------------------------------*//*
*//*[0]PETSC ERROR: Caught signal number 11 SEGV: Segmentation 
Violation, probably memory access out of range*//*
*//*[0]PETSC ERROR: Try option -start_in_debugger or 
-on_error_attach_debugger*//*
*//*[0]PETSC ERROR: or see 
http://www.mcs.anl.gov/petsc/documentation/faq.html#valgrind[0]PETSC 
ERROR: or try http://valgrind.org on GNU/linux and Apple Mac OS X to 
find memory corruption errors*//*
*//*[0]PETSC ERROR: likely location of problem given in stack below*//*
*//*[0]PETSC ERROR: ---------------------  Stack Frames 
------------------------------------*//*
*//*[0]PETSC ERROR: Note: The EXACT line numbers in the stack are not 
available,*//*
*//*[0]PETSC ERROR:       INSTEAD the line number of the start of the 
function*//*
*//*[0]PETSC ERROR:       is given.*//*
*//*[0]PETSC ERROR: --------------------- Error Message 
------------------------------------*//*
*//*[0]PETSC ERROR: Signal received!*//*
*//*[0]PETSC ERROR: 
------------------------------------------------------------------------*//*
*//*[0]PETSC ERROR: Petsc Release Version 3.4.4, Mar, 13, 2014 *//*
*//*[0]PETSC ERROR: See docs/changes/index.html for recent updates.*//*
*//*[0]PETSC ERROR: See docs/faq.html for hints about trouble shooting.*//*
*//*[0]PETSC ERROR: See docs/index.html for manual pages.*//*
*//*[0]PETSC ERROR: 
------------------------------------------------------------------------*//*
*//*[0]PETSC ERROR: ./*//*PgmPrincipal*//*on a arch-linux2-c-debug named 
r10dawesr by ndengues Fri Apr 25 09:36:59 2014*//*
*//*[0]PETSC ERROR: Libraries linked from 
/usr/local/home/ndengues/Downloads/Libraries/petsc-3.4.4/arch-linux2-c-debug/lib*//*
*//*[0]PETSC ERROR: Configure run at Thu Apr 24 20:35:39 2014*//*
*//*[0]PETSC ERROR: Configure options --with-cc=gcc --with-fc=gfortran 
--download-f-blas-lapack --download-mpich*//*
*//*[0]PETSC ERROR: 
------------------------------------------------------------------------*//*
*//*[0]PETSC ERROR: User provided function() line 0 in unknown directory 
unknown file*//*
*//*application called MPI_Abort(MPI_COMM_WORLD, 59) - process 0*//*
*//*[cli_0]: aborting job:*//*
*//*application called MPI_Abort(MPI_COMM_WORLD, 59) - process 0*/
***



Once more, the program execute without troubles for all the examples in 
the package.

Sincerely,


On 04/24/2014 10:14 PM, Satish Balay wrote:
> petsc examples work fine.
>
> send us the code that breaks. Also copy/paste the *complete* output from 'make' when
> you attempt to compile this code.
>
> Satish
>
> On Thu, 24 Apr 2014, Barry Smith wrote:
>
>>    Try .F files instead of .F90 some Fortran compilers handle preprocessing in different ways for .F and .F90
>>
>>     Please also send configure.log and make.log so we know what compilers etc you are using.
>>
>>    Barry
>>
>>    Note that we are using C pre processing directives,
>>
>> On Apr 24, 2014, at 7:57 PM, Steve Ndengue <steve.ndengue at gmail.com> wrote:
>>
>>> The results is the same with lower case letters...
>>>
>>>
>>> On 04/24/2014 07:39 PM, Matthew Knepley wrote:
>>>> On Thu, Apr 24, 2014 at 8:36 PM, Steve Ndengue <steve.ndengue at gmail.com> wrote:
>>>> Thank you for the quick reply.
>>>>
>>>> It partly solved the problem; the PETSC and SLEPC files are now included in the compilation.
>>>>
>>>> However the Warning are now error messages!
>>>> ***
>>>> pgm_hatom_offaxis_cyl.F90:880:0: error: invalid preprocessing directive #INCLUDE
>>>> pgm_hatom_offaxis_cyl.F90:887:0: error: invalid preprocessing directive #IF
>>>> pgm_hatom_offaxis_cyl.F90:890:0: error: invalid preprocessing directive #ELSE
>>>> pgm_hatom_offaxis_cyl.F90:893:0: error: invalid preprocessing directive #ENDIF
>>>> ***
>>>>
>>>> These should not be capitalized.
>>>>
>>>>    Matt
>>>>   
>>>> And the corresponding code lines are:
>>>> ***
>>>> #INCLUDE <finclude/slepcepsdef.h>
>>>>          USE slepceps
>>>>
>>>> !
>>>>   IMPLICIT NONE
>>>> !----------------------
>>>> !
>>>> #IF DEFINED(PETSC_USE_FORTRAN_DATATYPES)
>>>>          TYPE(Mat)                                          A
>>>>          TYPE(EPS)                                          solver
>>>> #ELSE
>>>>          Mat                                                A
>>>>          EPS                                                solver
>>>> #ENDIF
>>>> ***
>>>>
>>>> Sincerely.
>>>>
>>>> On 04/24/2014 07:25 PM, Satish Balay wrote:
>>>>> On Thu, 24 Apr 2014, Steve Ndengue wrote:
>>>>>
>>>>>
>>>>>> Dear all,
>>>>>>
>>>>>> I am having trouble compiling a code with some dependencies that shall call
>>>>>> SLEPC.
>>>>>> The compilation however goes perfectly for the various tutorials and tests in
>>>>>> the package.
>>>>>>
>>>>>> A sample makefile looks like:
>>>>>>
>>>>>> ***
>>>>>>
>>>>>> /default: code//
>>>>>> //routines: Rout1.o Rout2.o Rout3.o Module1.o//
>>>>>> //code: PgmPrincipal//
>>>>>> //sources=Rout1.f Rout2.f Rout3.f Module1.f90//
>>>>>>
>>>>> Its best to rename your sourcefiles .F/.F90 [and not .f/.f90
>>>>> [this enables using default targets from petsc makefiles - and its
>>>>> the correct notation for preprocessing - which is required by petsc/slepc]
>>>>>
>>>>>
>>>>>> //objets=Rout1.o Rout2.o Rout3.o Module1.o Pgm//Principal//.o//
>>>>>> //
>>>>>> //%.o: %.f//
>>>>>> //        -${FLINKER} -c $<//
>>>>>>
>>>>> And you shouldn't need to create the above .o.f target.
>>>>>
>>>>> Satish
>>>>>
>>>>>
>>>>>
>>>>>> //
>>>>>> //Module1_mod.mod Module1.o: //Module1//.f90//
>>>>>> //        -${FLINKER} -c //Module1//.f90//
>>>>>> //
>>>>>> //include ${SLEPC_DIR}/conf/slepc_common//
>>>>>> //
>>>>>> //Pgm//Principal//: ${objets} chkopts//
>>>>>> //        -${FLINKER} -o $@ ${objets} ${SLEPC_LIB}//
>>>>>> //
>>>>>> //.PHONY: clean//
>>>>>> //        ${RM} *.o *.mod Pgm//Principal//
>>>>>> /
>>>>>> ***
>>>>>>
>>>>>> The code exits with Warning and error messages of the type:
>>>>>>
>>>>>> *Warning: Pgm**Principal**.f90:880: Illegal preprocessor directive**
>>>>>> **Warning: Pgm**Principal**.f90:889: Illegal preprocessor directive**
>>>>>> **Warning: Pgm**Principal**.f90:892: Illegal preprocessor directive**
>>>>>> **Warning: Pgm**Principal**.f90:895: Illegal preprocessor directive*
>>>>>>
>>>>>> and
>>>>>>
>>>>>> *USE slepceps**
>>>>>> **                    1**
>>>>>> **Fatal Error: Can't open module file 'slepceps.mod' for reading at (1): No
>>>>>> such file or directory
>>>>>>
>>>>>>
>>>>>> *I do not get these errors with the tests and errors files.
>>>>>>
>>>>>>
>>>>>> Sincerely,
>>>>>>
>>>>>>
>>>>>>
>>>>
>>>> -- 
>>>> Steve
>>>>
>>>>
>>>>
>>>>
>>>> -- 
>>>> 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
>>>
>>> -- 
>>> Steve
>>>
>>>
>>>
>>


-- 
Steve

-------------- next part --------------
An HTML attachment was scrubbed...
URL: <http://lists.mcs.anl.gov/pipermail/petsc-users/attachments/20140425/200a1a0f/attachment.html>


More information about the petsc-users mailing list