Actual source code: zksp.c

  1: /*$Id: zksp.c,v 1.52 2001/08/07 21:32:16 bsmith Exp $*/

 3:  #include src/fortran/custom/zpetsc.h
 4:  #include petscksp.h

  6: #ifdef PETSC_HAVE_FORTRAN_CAPS
  7: #define kspgetresidualnorm_        KSPGETRESIDUALNORM
  8: #define kspgetconvergedreason_     KSPGETCONVERGEDREASON
  9: #define kspfgmressetmodifypc_      KSPFGMRESSETMODIFYPC
 10: #define kspfgmresmodifypcsles_     KSPFGMRESMODIFYPCSLES
 11: #define kspfgmresmodifypcnochange_ KSPFGMRESMODIFYPCNOCHANGE
 12: #define kspdefaultconverged_       KSPDEFAULTCONVERGED
 13: #define kspskipconverged_          KSPSKIPCONVERGED
 14: #define kspgmreskrylovmonitor_     KSPGMRESKRYLOVMONITOR
 15: #define kspdefaultmonitor_         KSPDEFAULTMONITOR
 16: #define ksptruemonitor_            KSPTRUEMONITOR
 17: #define kspvecviewmonitor_         KSPVECVIEWMONITOR
 18: #define ksplgmonitor_              KSPLGMONITOR
 19: #define ksplgtruemonitor_          KSPLGTRUEMONITOR
 20: #define kspsingularvaluemonitor_   KSPSINGULARVALUEMONITOR
 21: #define kspregisterdestroy_        KSPREGISTERDESTROY
 22: #define kspdestroy_                KSPDESTROY
 23: #define ksplgmonitordestroy_       KSPLGMONITORDESTROY
 24: #define ksplgmonitorcreate_        KSPLGMONITORCREATE
 25: #define kspgetrhs_                 KSPGETRHS
 26: #define kspgetsolution_            KSPGETSOLUTION
 27: #define kspgetpc_                  KSPGETPC
 28: #define kspsetmonitor_             KSPSETMONITOR
 29: #define kspsetconvergencetest_     KSPSETCONVERGENCETEST
 30: #define kspcreate_                 KSPCREATE
 31: #define kspsetoptionsprefix_       KSPSETOPTIONSPREFIX
 32: #define kspappendoptionsprefix_    KSPAPPENDOPTIONSPREFIX
 33: #define kspgettype_                KSPGETTYPE
 34: #define kspgetpreconditionerside_  KSPGETPRECONDITIONERSIDE
 35: #define kspbuildsolution_          KSPBUILDSOLUTION 
 36: #define kspbuildresidual_          KSPBUILDRESIDUAL
 37: #define kspsettype_                KSPSETTYPE           
 38: #define kspgetresidualhistory_     KSPGETRESIDUALHISTORY
 39: #define kspgetoptionsprefix_       KSPGETOPTIONSPREFIX
 40: #define kspview_                   KSPVIEW
 41: #define kspgmressetrestart_        KSPGMRESSETRESTART
 42: #define kspsetnormtype_            KSPSETNORMTYPE
 43: #elif !defined(PETSC_HAVE_FORTRAN_UNDERSCORE)
 44: #define kspgetconvergedreason_     kspgetconvergedreason
 45: #define kspfgmressetmodifypc_      kspfgmressetmodifypc
 46: #define kspfgmresmodifypcsles_     kspfgmresmodifypcsles
 47: #define kspfgmresmodifypcnochange_ kspfgmresmodifypcnochange
 48: #define kspdefaultconverged_       kspdefaultconverged
 49: #define kspskipconverged_          kspskipconverged
 50: #define kspsingularvaluemonitor_   kspsingularvaluemonitor
 51: #define kspgmreskrylovmonitor_     kspgmreskrylovmonitor
 52: #define kspdefaultmonitor_         kspdefaultmonitor
 53: #define ksptruemonitor_            ksptruemonitor
 54: #define kspvecviewmonitor_         kspvecviewmonitor
 55: #define ksplgmonitor_              ksplgmonitor
 56: #define ksplgtruemonitor_          ksplgtruemonitor
 57: #define kspgetresidualhistory_     kspgetresidualhistory
 58: #define kspsettype_                kspsettype
 59: #define kspregisterdestroy_        kspregisterdestroy
 60: #define kspdestroy_                kspdestroy
 61: #define ksplgmonitordestroy_       ksplgmonitordestroy
 62: #define ksplgmonitorcreate_        ksplgmonitorcreate
 63: #define kspgetrhs_                 kspgetrhs
 64: #define kspgetsolution_            kspgetsolution
 65: #define kspgetpc_                  kspgetpc
 66: #define kspsetmonitor_             kspsetmonitor
 67: #define kspsetconvergencetest_     kspsetconvergencetest
 68: #define kspcreate_                 kspcreate
 69: #define kspsetoptionsprefix_       kspsetoptionsprefix
 70: #define kspappendoptionsprefix_    kspappendoptionsprefix
 71: #define kspgettype_                kspgettype
 72: #define kspgetpreconditionerside_  kspgetpreconditionerside
 73: #define kspbuildsolution_          kspbuildsolution
 74: #define kspbuildresidual_          kspbuildresidual
 75: #define kspgetoptionsprefix_       kspgetoptionsprefix
 76: #define kspview_                   kspview
 77: #define kspgetresidualnorm_        kspgetresidualnorm
 78: #define kspgmressetrestart_        kspgmressetrestart
 79: #define kspsetnormtype_            kspsetnormtype
 80: #endif

 82: EXTERN_C_BEGIN

 84: void PETSC_STDCALL kspgmressetrestart_(KSP *ksp,int *max_k, int *ierr )
 85: {
 86:   *KSPGMRESSetRestart(*ksp,*max_k);
 87: }

 89: void PETSC_STDCALL kspgetresidualnorm_(KSP *ksp,PetscReal *rnorm,int *ierr)
 90: {
 91:   *KSPGetResidualNorm(*ksp,rnorm);
 92: }

 94: void PETSC_STDCALL kspgetconvergedreason_(KSP *ksp,KSPConvergedReason *reason,int *ierr)
 95: {
 96:   *KSPGetConvergedReason(*ksp,reason);
 97: }

 99: /* function */
100: void PETSC_STDCALL kspview_(KSP *ksp,PetscViewer *viewer, int *ierr)
101: {
102:   PetscViewer v;
103:   PetscPatchDefaultViewers_Fortran(viewer,v);
104:   *KSPView(*ksp,v);
105: }

107: void kspdefaultconverged_(KSP *ksp,int *n,PetscReal *rnorm,KSPConvergedReason *flag,void *dummy,int *ierr)
108: {
109:   CHKFORTRANNULLOBJECT(dummy);
110:   *KSPDefaultConverged(*ksp,*n,*rnorm,flag,dummy);
111: }

113: void kspskipconverged_(KSP *ksp,int *n,PetscReal *rnorm,KSPConvergedReason *flag,void *dummy,int *ierr)
114: {
115:   CHKFORTRANNULLOBJECT(dummy);
116:   *KSPSkipConverged(*ksp,*n,*rnorm,flag,dummy);
117: }

119: void PETSC_STDCALL kspgetresidualhistory_(KSP *ksp,int *na,int *ierr)
120: {
121:   *KSPGetResidualHistory(*ksp,PETSC_NULL,na);
122: }

124: void PETSC_STDCALL kspsettype_(KSP *ksp,CHAR type PETSC_MIXED_LEN(len),int *ierr PETSC_END_LEN(len))
125: {
126:   char *t;

128:   FIXCHAR(type,len,t);
129:   *KSPSetType(*ksp,t);
130:   FREECHAR(type,t);
131: }

133: void PETSC_STDCALL kspgettype_(KSP *ksp,CHAR name PETSC_MIXED_LEN(len),int *ierr PETSC_END_LEN(len))
134: {
135:   char *tname;

137:   *KSPGetType(*ksp,&tname);if (*ierr) return;
138: #if defined(PETSC_USES_CPTOFCD)
139:   {
140:     char *t = _fcdtocp(name); int len1 = _fcdlen(name);
141:     *PetscStrncpy(t,tname,len1);
142:   }
143: #else
144:   *PetscStrncpy(name,tname,len);
145: #endif
146: }

148: void PETSC_STDCALL kspgetpreconditionerside_(KSP *ksp,PCSide *side,int *ierr){
149: *KSPGetPreconditionerSide(*ksp,side);
150: }

152: void PETSC_STDCALL kspsetoptionsprefix_(KSP *ksp,CHAR prefix PETSC_MIXED_LEN(len),
153:                                         int *ierr PETSC_END_LEN(len))
154: {
155:   char *t;

157:   FIXCHAR(prefix,len,t);
158:   *KSPSetOptionsPrefix(*ksp,t);
159:   FREECHAR(prefix,t);
160: }

162: void PETSC_STDCALL kspappendoptionsprefix_(KSP *ksp,CHAR prefix PETSC_MIXED_LEN(len),
163:                                            int *ierr PETSC_END_LEN(len))
164: {
165:   char *t;

167:   FIXCHAR(prefix,len,t);
168:   *KSPAppendOptionsPrefix(*ksp,t);
169:   FREECHAR(prefix,t);
170: }

172: void PETSC_STDCALL kspcreate_(MPI_Comm *comm,KSP *ksp,int *ierr){
173:   *KSPCreate((MPI_Comm)PetscToPointerComm(*comm),ksp);
174: }

176: static void (PETSC_STDCALL *f2)(KSP*,int*,PetscReal*,KSPConvergedReason*,void*,int*);
177: static int ourtest(KSP ksp,int i,PetscReal d,KSPConvergedReason *reason,void* ctx)
178: {
180:   (*f2)(&ksp,&i,&d,reason,ctx,&ierr);
181:   return 0;
182: }
183: void PETSC_STDCALL kspsetconvergencetest_(KSP *ksp,
184:       void (PETSC_STDCALL *converge)(KSP*,int*,PetscReal*,KSPConvergedReason*,void*,int*),void *cctx,int *ierr)
185: {
186:   if ((void(*)(void))converge == (void(*)(void))kspdefaultconverged_) {
187:     *KSPSetConvergenceTest(*ksp,KSPDefaultConverged,0);
188:   } else if ((void(*)(void))converge == (void(*)(void))kspskipconverged_) {
189:     *KSPSetConvergenceTest(*ksp,KSPSkipConverged,0);
190:   } else {
191:     f2 = converge;
192:     *KSPSetConvergenceTest(*ksp,ourtest,cctx);
193:   }
194: }

196: /*
197:         These are not usually called from Fortran but allow Fortran users 
198:    to transparently set these monitors from .F code
199:    
200:    functions, hence no STDCALL
201: */
202: void kspgmreskrylovmonitor_(KSP *ksp,int *it,PetscReal *norm,void *ctx,int *ierr)
203: {
204:   *KSPGMRESKrylovMonitor(*ksp,*it,*norm,ctx);
205: }

207: void  kspdefaultmonitor_(KSP *ksp,int *it,PetscReal *norm,void *ctx,int *ierr)
208: {
209:   *KSPDefaultMonitor(*ksp,*it,*norm,ctx);
210: }
211: 
212: void  kspsingularvaluemonitor_(KSP *ksp,int *it,PetscReal *norm,void *ctx,int *ierr)
213: {
214:   *KSPSingularValueMonitor(*ksp,*it,*norm,ctx);
215: }

217: void  ksplgmonitor_(KSP *ksp,int *it,PetscReal *norm,void *ctx,int *ierr)
218: {
219:   *KSPLGMonitor(*ksp,*it,*norm,ctx);
220: }

222: void  ksplgtruemonitor_(KSP *ksp,int *it,PetscReal *norm,void *ctx,int *ierr)
223: {
224:   *KSPLGTrueMonitor(*ksp,*it,*norm,ctx);
225: }

227: void  ksptruemonitor_(KSP *ksp,int *it,PetscReal *norm,void *ctx,int *ierr)
228: {
229:   *KSPTrueMonitor(*ksp,*it,*norm,ctx);
230: }

232: void  kspvecviewmonitor_(KSP *ksp,int *it,PetscReal *norm,void *ctx,int *ierr)
233: {
234:   *KSPVecViewMonitor(*ksp,*it,*norm,ctx);
235: }

237: static void (PETSC_STDCALL *f1)(KSP*,int*,PetscReal*,void*,int*);
238: static int ourmonitor(KSP ksp,int i,PetscReal d,void* ctx)
239: {
240:   int 0;
241:   (*f1)(&ksp,&i,&d,ctx,&ierr);
242:   return 0;
243: }
244: static void (PETSC_STDCALL *f21)(void*,int*);
245: static int ourdestroy(void* ctx)
246: {
247:   int 0;
248:   (*f21)(ctx,&ierr);
249:   return 0;
250: }

252: void PETSC_STDCALL kspsetmonitor_(KSP *ksp,void (PETSC_STDCALL *monitor)(KSP*,int*,PetscReal*,void*,int*),
253:                     void *mctx,void (PETSC_STDCALL *monitordestroy)(void *,int *),int *ierr)
254: {
255:   if ((void(*)(void))monitor == (void(*)(void))kspdefaultmonitor_) {
256:     *KSPSetMonitor(*ksp,KSPDefaultMonitor,0,0);
257:   } else if ((void(*)(void))monitor == (void(*)(void))ksplgmonitor_) {
258:     *KSPSetMonitor(*ksp,KSPLGMonitor,0,0);
259:   } else if ((void(*)(void))monitor == (void(*)(void))ksplgtruemonitor_) {
260:     *KSPSetMonitor(*ksp,KSPLGTrueMonitor,0,0);
261:   } else if ((void(*)(void))monitor == (void(*)(void))kspvecviewmonitor_) {
262:     *KSPSetMonitor(*ksp,KSPVecViewMonitor,0,0);
263:   } else if ((void(*)(void))monitor == (void(*)(void))ksptruemonitor_) {
264:     *KSPSetMonitor(*ksp,KSPTrueMonitor,0,0);
265:   } else if ((void(*)(void))monitor == (void(*)(void))kspsingularvaluemonitor_) {
266:     *KSPSetMonitor(*ksp,KSPSingularValueMonitor,0,0);
267:   } else {
268:     f1  = monitor;
269:     if (FORTRANNULLFUNCTION(monitordestroy)) {
270:       *KSPSetMonitor(*ksp,ourmonitor,mctx,0);
271:     } else {
272:       f21 = monitordestroy;
273:       *KSPSetMonitor(*ksp,ourmonitor,mctx,ourdestroy);
274:     }
275:   }
276: }

278: void PETSC_STDCALL kspgetpc_(KSP *ksp,PC *B,int *ierr)
279: {
280:   *KSPGetPC(*ksp,B);
281: }

283: void PETSC_STDCALL kspgetsolution_(KSP *ksp,Vec *v,int *ierr)
284: {
285:   *KSPGetSolution(*ksp,v);
286: }

288: void PETSC_STDCALL kspgetrhs_(KSP *ksp,Vec *r,int *ierr)
289: {
290:   *KSPGetRhs(*ksp,r);
291: }

293: /*
294:    Possible bleeds memory but cannot be helped.
295: */
296: void PETSC_STDCALL ksplgmonitorcreate_(CHAR host PETSC_MIXED_LEN(len1),
297:                     CHAR label PETSC_MIXED_LEN(len2),int *x,int *y,int *m,int *n,PetscDrawLG *ctx,
298:                     int *ierr PETSC_END_LEN(len1) PETSC_END_LEN(len2))
299: {
300:   char   *t1,*t2;

302:   FIXCHAR(host,len1,t1);
303:   FIXCHAR(label,len2,t2);
304:   *KSPLGMonitorCreate(t1,t2,*x,*y,*m,*n,ctx);
305: }

307: void PETSC_STDCALL ksplgmonitordestroy_(PetscDrawLG *ctx,int *ierr)
308: {
309:   *KSPLGMonitorDestroy(*ctx);
310: }

312: void PETSC_STDCALL kspdestroy_(KSP *ksp,int *ierr)
313: {
314:   *KSPDestroy(*ksp);
315: }

317: void PETSC_STDCALL kspregisterdestroy_(int* ierr)
318: {
319:   *KSPRegisterDestroy();
320: }

322: void PETSC_STDCALL kspbuildsolution_(KSP *ctx,Vec *v,Vec *V,int *ierr)
323: {
324:   *KSPBuildSolution(*ctx,*v,V);
325: }

327: void PETSC_STDCALL kspbuildresidual_(KSP *ctx,Vec *t,Vec *v,Vec *V,int *ierr)
328: {
329:   *KSPBuildResidual(*ctx,*t,*v,V);
330: }

332: void PETSC_STDCALL kspgetoptionsprefix_(KSP *ksp,CHAR prefix PETSC_MIXED_LEN(len),int *ierr PETSC_END_LEN(len))
333: {
334:   char *tname;

336:   *KSPGetOptionsPrefix(*ksp,&tname);
337: #if defined(PETSC_USES_CPTOFCD)
338:   {
339:     char *t = _fcdtocp(prefix); int len1 = _fcdlen(prefix);
340:     *PetscStrncpy(t,tname,len1); if (*ierr) return;
341:   }
342: #else
343:   *PetscStrncpy(prefix,tname,len); if (*ierr) return;
344: #endif
345: }

347: static void (PETSC_STDCALL *f109)(KSP*,int*,int*,PetscReal*,void*,int*);
348: static int ourmodify(KSP ksp,int i,int i2,PetscReal d,void* ctx)
349: {
350:   int 0;
351:   (*f109)(&ksp,&i,&i2,&d,ctx,&ierr);
352:   return 0;
353: }

355: static void (PETSC_STDCALL *f210)(void*,int*);
356: static int ourmoddestroy(void* ctx)
357: {
358:   int 0;
359:   (*f210)(ctx,&ierr);
360:   return 0;
361: }

363: void PETSC_STDCALL kspfgmresmodifypcnochange_(KSP *ksp,int *total_its,int *loc_its,PetscReal *res_norm,void* dummy,int *ierr)
364: {
365:   *KSPFGMRESModifyPCNoChange(*ksp,*total_its,*loc_its,*res_norm,dummy);
366: }

368: void PETSC_STDCALL kspfgmresmodifypcsles_(KSP *ksp,int *total_its,int *loc_its,PetscReal *res_norm,void*dummy,int *ierr)
369: {
370:   *KSPFGMRESModifyPCSLES(*ksp,*total_its,*loc_its,*res_norm,dummy);
371: }

373: void PETSC_STDCALL kspfgmressetmodifypc_(KSP *ksp,void (PETSC_STDCALL *fcn)(KSP*,int*,int*,PetscReal*,void*,int*),void* ctx,void (PETSC_STDCALL *d)(void*,int*),int *ierr)
374: {
375:   if ((void(*)(void))fcn == (void(*)(void))kspfgmresmodifypcsles_) {
376:     *KSPFGMRESSetModifyPC(*ksp,KSPFGMRESModifyPCSLES,0,0);
377:   } else if ((void(*)(void))fcn == (void(*)(void))kspfgmresmodifypcnochange_) {
378:     *KSPFGMRESSetModifyPC(*ksp,KSPFGMRESModifyPCNoChange,0,0);
379:   } else {
380:     f109 = fcn;
381:     if (FORTRANNULLFUNCTION(d)) {
382:       *KSPFGMRESSetModifyPC(*ksp,ourmodify,ctx,0);
383:     } else {
384:       f210 = d;
385:       *KSPFGMRESSetModifyPC(*ksp,ourmodify,ctx,ourmoddestroy);
386:     }
387:   }
388: }

390: void PETSC_STDCALL kspsetnormtype_(KSP *ksp,KSPNormType *type,int *ierr)
391: {
392:   *KSPSetNormType(*ksp,*type);
393: }

395: EXTERN_C_END