xref: /petsc/src/sys/error/pstack.c (revision 16ad03002c971217bb7d68d241ce25e1b2bc6bef)
1d9262e54SJed Brown 
2c6db04a5SJed Brown #include <petscsys.h>        /*I  "petscsys.h"   I*/
3d9262e54SJed Brown 
48bf1f09cSShri Abhyankar #if defined(PETSC_USE_DEBUG)
5d9262e54SJed Brown 
61f46d60fSShri Abhyankar #if defined(PETSC_HAVE_PTHREADCLASSES)
7997ce2baSShri Abhyankar #if defined(PETSC_PTHREAD_LOCAL)
8997ce2baSShri Abhyankar PETSC_PTHREAD_LOCAL PetscStack *petscstack = 0;
91f46d60fSShri Abhyankar #endif
10fe89fe5aSShri Abhyankar #else
117087cfbeSBarry Smith PetscStack *petscstack = 0;
12fe89fe5aSShri Abhyankar #endif
13d9262e54SJed Brown 
1415681b3cSBarry Smith 
15e04113cfSBarry Smith #if defined(PETSC_HAVE_SAWS)
16e04113cfSBarry Smith #include <petscviewersaws.h>
1715681b3cSBarry Smith 
182657e9d9SBarry Smith static PetscBool amsmemstack = PETSC_FALSE;
1915681b3cSBarry Smith 
20d9262e54SJed Brown #undef __FUNCT__
21e04113cfSBarry Smith #define __FUNCT__ "PetscStackSAWsGrantAccess"
2215681b3cSBarry Smith /*@C
23e04113cfSBarry Smith    PetscStackSAWsGrantAccess - Grants access of the PETSc stack frames to the SAWs publisher
2415681b3cSBarry Smith 
2515681b3cSBarry Smith    Collective on PETSC_COMM_WORLD?
2615681b3cSBarry Smith 
2715681b3cSBarry Smith    Level: developer
2815681b3cSBarry Smith 
2915681b3cSBarry Smith    Concepts: publishing object
3015681b3cSBarry Smith 
31e04113cfSBarry Smith    Developers Note: Cannot use PetscFunctionBegin/Return() or PetscStackCallSAWs() since it may be used within those routines
3215681b3cSBarry Smith 
33e04113cfSBarry Smith .seealso: PetscObjectSetName(), PetscObjectSAWsViewOff(), PetscObjectSAWsTakeAccess()
3415681b3cSBarry Smith 
3515681b3cSBarry Smith @*/
36e04113cfSBarry Smith void  PetscStackSAWsGrantAccess(void)
37d9262e54SJed Brown {
38ec957eceSBarry Smith   if (amsmemstack) {
39*16ad0300SBarry Smith     /* ignore any errors from SAWs */
409a492a5cSBarry Smith     SAWs_Unlock();
4115681b3cSBarry Smith   }
42d9262e54SJed Brown }
43d9262e54SJed Brown 
44d9262e54SJed Brown #undef __FUNCT__
45e04113cfSBarry Smith #define __FUNCT__ "PetscStackSAWsTakeAccess"
4615681b3cSBarry Smith /*@C
47e04113cfSBarry Smith    PetscStackSAWsTakeAccess - Takes access of the PETSc stack frames to the SAWs publisher
4815681b3cSBarry Smith 
4915681b3cSBarry Smith    Collective on PETSC_COMM_WORLD?
5015681b3cSBarry Smith 
5115681b3cSBarry Smith    Level: developer
5215681b3cSBarry Smith 
5315681b3cSBarry Smith    Concepts: publishing object
5415681b3cSBarry Smith 
55e04113cfSBarry Smith    Developers Note: Cannot use PetscFunctionBegin/Return() or PetscStackCallSAWs() since it may be used within those routines
5615681b3cSBarry Smith 
57e04113cfSBarry Smith .seealso: PetscObjectSetName(), PetscObjectSAWsViewOff(), PetscObjectSAWsTakeAccess()
5815681b3cSBarry Smith 
5915681b3cSBarry Smith @*/
60e04113cfSBarry Smith void  PetscStackSAWsTakeAccess(void)
61d9262e54SJed Brown {
62ec957eceSBarry Smith   if (amsmemstack) {
63*16ad0300SBarry Smith     /* ignore any errors from SAWs */
649a492a5cSBarry Smith     SAWs_Lock();
6515681b3cSBarry Smith   }
6615681b3cSBarry Smith }
6715681b3cSBarry Smith 
68e04113cfSBarry Smith PetscErrorCode PetscStackViewSAWs(void)
6915681b3cSBarry Smith {
7015681b3cSBarry Smith   PetscStack*    petscstackp;
7115681b3cSBarry Smith 
7215681b3cSBarry Smith   petscstackp = (PetscStack*)PetscThreadLocalGetValue(petscstack);
7367a1e3b6SBarry Smith   PetscStackCallSAWs(SAWs_Register,("/PETSc/Stack/functions",petscstackp->function,20,SAWs_READ,SAWs_STRING));
74b8a862caSBarry Smith   PetscStackCallSAWs(SAWs_Register,("/PETSc/Stack/_current_size",&petscstackp->currentsize,1,SAWs_READ,SAWs_INT));
752657e9d9SBarry Smith   amsmemstack = PETSC_TRUE;
7615681b3cSBarry Smith   return 0;
7715681b3cSBarry Smith }
7815681b3cSBarry Smith 
7915681b3cSBarry Smith #undef __FUNCT__
80e04113cfSBarry Smith #define __FUNCT__ "PetscStackSAWsViewOff"
81e04113cfSBarry Smith PetscErrorCode PetscStackSAWsViewOff(void)
8215681b3cSBarry Smith {
83d9262e54SJed Brown   PetscFunctionBegin;
84*16ad0300SBarry Smith   PetscStackCallSAWs(SAWs_Delete,("/PETSc/Stack"));
852657e9d9SBarry Smith   amsmemstack = PETSC_FALSE;
86d9262e54SJed Brown   PetscFunctionReturn(0);
87d9262e54SJed Brown }
88d9262e54SJed Brown 
8915681b3cSBarry Smith #endif
9015681b3cSBarry Smith 
917d5f7e0cSShri Abhyankar PetscErrorCode PetscStackCreate(void)
9274b43855SShri Abhyankar {
9374b43855SShri Abhyankar   PetscStack *petscstack_in;
94dbf62e16SBarry Smith   if (PetscStackActive()) return 0;
9574b43855SShri Abhyankar 
9674b43855SShri Abhyankar   petscstack_in              = (PetscStack*)malloc(sizeof(PetscStack));
9774b43855SShri Abhyankar   petscstack_in->currentsize = 0;
98047240e1SBarry Smith   PetscThreadLocalSetValue((PetscThreadKey*)&petscstack,petscstack_in);
9915681b3cSBarry Smith 
100e04113cfSBarry Smith #if defined(PETSC_HAVE_SAWS)
10115681b3cSBarry Smith   {
10215681b3cSBarry Smith   PetscBool flg = PETSC_FALSE;
10315681b3cSBarry Smith   PetscOptionsHasName(NULL,"-stack_view",&flg);
104e04113cfSBarry Smith   if (flg) PetscStackViewSAWs();
10515681b3cSBarry Smith   }
10615681b3cSBarry Smith #endif
10774b43855SShri Abhyankar   return 0;
10874b43855SShri Abhyankar }
10974b43855SShri Abhyankar 
11015681b3cSBarry Smith 
111d9262e54SJed Brown #undef __FUNCT__
112d9262e54SJed Brown #define __FUNCT__ "PetscStackView"
113639ff905SBarry Smith PetscErrorCode  PetscStackView(FILE *file)
114d9262e54SJed Brown {
115d9262e54SJed Brown   int        i;
1161f46d60fSShri Abhyankar   PetscStack *petscstackp;
117d9262e54SJed Brown 
1181f46d60fSShri Abhyankar   petscstackp = (PetscStack*)PetscThreadLocalGetValue(petscstack);
119639ff905SBarry Smith   if (!file) file = PETSC_STDOUT;
120d9262e54SJed Brown 
121d9262e54SJed Brown   if (file == PETSC_STDOUT) {
122d9262e54SJed Brown     (*PetscErrorPrintf)("Note: The EXACT line numbers in the stack are not available,\n");
123d9262e54SJed Brown     (*PetscErrorPrintf)("      INSTEAD the line number of the start of the function\n");
124d9262e54SJed Brown     (*PetscErrorPrintf)("      is given.\n");
125a297a907SKarl Rupp     for (i=petscstackp->currentsize-1; i>=0; i--) (*PetscErrorPrintf)("[%d] %s line %d %s%s\n",PetscGlobalRank,petscstackp->function[i],petscstackp->line[i],petscstackp->directory[i],petscstackp->file[i]);
126d9262e54SJed Brown   } else {
127d9262e54SJed Brown     fprintf(file,"Note: The EXACT line numbers in the stack are not available,\n");
128d9262e54SJed Brown     fprintf(file,"      INSTEAD the line number of the start of the function\n");
129d9262e54SJed Brown     fprintf(file,"      is given.\n");
130a297a907SKarl Rupp     for (i=petscstackp->currentsize-1; i>=0; i--) fprintf(file,"[%d] %s line %d %s%s\n",PetscGlobalRank,petscstackp->function[i],petscstackp->line[i],petscstackp->directory[i],petscstackp->file[i]);
131d9262e54SJed Brown   }
132d9262e54SJed Brown   return 0;
133d9262e54SJed Brown }
134d9262e54SJed Brown 
1357d5f7e0cSShri Abhyankar PetscErrorCode PetscStackDestroy(void)
13674b43855SShri Abhyankar {
137dbf62e16SBarry Smith   if (PetscStackActive()) {
1381f46d60fSShri Abhyankar     PetscStack *petscstack_in;
1391f46d60fSShri Abhyankar     petscstack_in = (PetscStack*)PetscThreadLocalGetValue(petscstack);
14074b43855SShri Abhyankar     free(petscstack_in);
141047240e1SBarry Smith     PetscThreadLocalSetValue((PetscThreadKey*)&petscstack,(PetscStack*)0);
14261d886c9SShri Abhyankar     PetscThreadLocalDestroy(petscstack); /* Deletes pthread_key if it was used */
1437d5f7e0cSShri Abhyankar   }
144d9262e54SJed Brown   return 0;
145d9262e54SJed Brown }
146d9262e54SJed Brown 
147d9262e54SJed Brown #undef __FUNCT__
148d9262e54SJed Brown #define __FUNCT__ "PetscStackCopy"
149d9262e54SJed Brown /*  PetscFunctionBegin;  so that make rule checkbadPetscFunctionBegin works */
1507087cfbeSBarry Smith PetscErrorCode  PetscStackCopy(PetscStack *sint,PetscStack *sout)
151d9262e54SJed Brown {
152d9262e54SJed Brown   int i;
153d9262e54SJed Brown 
154a297a907SKarl Rupp   if (!sint) sout->currentsize = 0;
155a297a907SKarl Rupp   else {
156d9262e54SJed Brown     for (i=0; i<sint->currentsize; i++) {
157d9262e54SJed Brown       sout->function[i]     = sint->function[i];
158d9262e54SJed Brown       sout->file[i]         = sint->file[i];
159d9262e54SJed Brown       sout->directory[i]    = sint->directory[i];
160d9262e54SJed Brown       sout->line[i]         = sint->line[i];
161a8d2bbe5SBarry Smith       sout->petscroutine[i] = sint->petscroutine[i];
162d9262e54SJed Brown     }
163d9262e54SJed Brown     sout->currentsize = sint->currentsize;
164d9262e54SJed Brown   }
165d9262e54SJed Brown   return 0;
166d9262e54SJed Brown }
167d9262e54SJed Brown 
168d9262e54SJed Brown #undef __FUNCT__
169d9262e54SJed Brown #define __FUNCT__ "PetscStackPrint"
170d9262e54SJed Brown /*  PetscFunctionBegin;  so that make rule checkbadPetscFunctionBegin works */
1717087cfbeSBarry Smith PetscErrorCode  PetscStackPrint(PetscStack *sint,FILE *fp)
172d9262e54SJed Brown {
173d9262e54SJed Brown   int i;
174d9262e54SJed Brown 
175d9262e54SJed Brown   if (!sint) return(0);
176a297a907SKarl Rupp   for (i=sint->currentsize-2; i>=0; i--) fprintf(fp,"      [%d]  %s() line %d in %s%s\n",PetscGlobalRank,sint->function[i],sint->line[i],sint->directory[i],sint->file[i]);
177d9262e54SJed Brown   return 0;
178d9262e54SJed Brown }
179d9262e54SJed Brown 
180d9262e54SJed Brown #else
18115681b3cSBarry Smith 
182d9262e54SJed Brown #undef __FUNCT__
183d9262e54SJed Brown #define __FUNCT__ "PetscStackCreate"
1847087cfbeSBarry Smith PetscErrorCode  PetscStackCreate(void)
185d9262e54SJed Brown {
186d9262e54SJed Brown   PetscFunctionBegin;
187d9262e54SJed Brown   PetscFunctionReturn(0);
188d9262e54SJed Brown }
189d9262e54SJed Brown #undef __FUNCT__
190d9262e54SJed Brown #define __FUNCT__ "PetscStackView"
191639ff905SBarry Smith PetscErrorCode  PetscStackView(FILE *file)
192d9262e54SJed Brown {
193d9262e54SJed Brown   PetscFunctionBegin;
194d9262e54SJed Brown   PetscFunctionReturn(0);
195d9262e54SJed Brown }
196d9262e54SJed Brown #undef __FUNCT__
197d9262e54SJed Brown #define __FUNCT__ "PetscStackDestroy"
1987087cfbeSBarry Smith PetscErrorCode  PetscStackDestroy(void)
199d9262e54SJed Brown {
200d9262e54SJed Brown   PetscFunctionBegin;
201d9262e54SJed Brown   PetscFunctionReturn(0);
202d9262e54SJed Brown }
203d9262e54SJed Brown 
204e04113cfSBarry Smith #if defined(PETSC_HAVE_SAWS)     /* SAWs stack functions do nothing in optimized mode */
205e04113cfSBarry Smith void PetscStackSAWsGrantAccess(void) {}
206e04113cfSBarry Smith void PetscStackSAWsTakeAccess(void) {}
2076eb7abf7SJed Brown 
208e04113cfSBarry Smith PetscErrorCode PetscStackViewSAWs(void)
2096eb7abf7SJed Brown {
2106eb7abf7SJed Brown   return 0;
2116eb7abf7SJed Brown }
2126eb7abf7SJed Brown 
2136eb7abf7SJed Brown #undef __FUNCT__
214e04113cfSBarry Smith #define __FUNCT__ "PetscStackSAWsViewOff"
215e04113cfSBarry Smith PetscErrorCode  PetscStackSAWsViewOff(void)
2166eb7abf7SJed Brown {
2176eb7abf7SJed Brown   PetscFunctionBegin;
2186eb7abf7SJed Brown   PetscFunctionReturn(0);
2196eb7abf7SJed Brown }
2206eb7abf7SJed Brown #endif
2216eb7abf7SJed Brown 
222d9262e54SJed Brown #endif
223d9262e54SJed Brown 
224