Actual source code: fnbasic.c

slepc-3.7.3 2016-09-29
Report Typos and Errors
  1: /*
  2:      Basic routines

  4:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  5:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  6:    Copyright (c) 2002-2016, Universitat Politecnica de Valencia, Spain

  8:    This file is part of SLEPc.

 10:    SLEPc is free software: you can redistribute it and/or modify it under  the
 11:    terms of version 3 of the GNU Lesser General Public License as published by
 12:    the Free Software Foundation.

 14:    SLEPc  is  distributed in the hope that it will be useful, but WITHOUT  ANY
 15:    WARRANTY;  without even the implied warranty of MERCHANTABILITY or  FITNESS
 16:    FOR  A  PARTICULAR PURPOSE. See the GNU Lesser General Public  License  for
 17:    more details.

 19:    You  should have received a copy of the GNU Lesser General  Public  License
 20:    along with SLEPc. If not, see <http://www.gnu.org/licenses/>.
 21:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 22: */

 24: #include <slepc/private/fnimpl.h>      /*I "slepcfn.h" I*/
 25: #include <slepcblaslapack.h>

 27: PetscFunctionList FNList = 0;
 28: PetscBool         FNRegisterAllCalled = PETSC_FALSE;
 29: PetscClassId      FN_CLASSID = 0;
 30: PetscLogEvent     FN_Evaluate = 0;
 31: static PetscBool  FNPackageInitialized = PETSC_FALSE;

 35: /*@C
 36:    FNFinalizePackage - This function destroys everything in the Slepc interface
 37:    to the FN package. It is called from SlepcFinalize().

 39:    Level: developer

 41: .seealso: SlepcFinalize()
 42: @*/
 43: PetscErrorCode FNFinalizePackage(void)
 44: {

 48:   PetscFunctionListDestroy(&FNList);
 49:   FNPackageInitialized = PETSC_FALSE;
 50:   FNRegisterAllCalled  = PETSC_FALSE;
 51:   return(0);
 52: }

 56: /*@C
 57:   FNInitializePackage - This function initializes everything in the FN package.
 58:   It is called from PetscDLLibraryRegister() when using dynamic libraries, and
 59:   on the first call to FNCreate() when using static libraries.

 61:   Level: developer

 63: .seealso: SlepcInitialize()
 64: @*/
 65: PetscErrorCode FNInitializePackage(void)
 66: {
 67:   char             logList[256];
 68:   char             *className;
 69:   PetscBool        opt;
 70:   PetscErrorCode   ierr;

 73:   if (FNPackageInitialized) return(0);
 74:   FNPackageInitialized = PETSC_TRUE;
 75:   /* Register Classes */
 76:   PetscClassIdRegister("Math Function",&FN_CLASSID);
 77:   /* Register Constructors */
 78:   FNRegisterAll();
 79:   /* Register Events */
 80:   PetscLogEventRegister("FNEvaluate",FN_CLASSID,&FN_Evaluate);
 81:   /* Process info exclusions */
 82:   PetscOptionsGetString(NULL,NULL,"-info_exclude",logList,256,&opt);
 83:   if (opt) {
 84:     PetscStrstr(logList,"fn",&className);
 85:     if (className) {
 86:       PetscInfoDeactivateClass(FN_CLASSID);
 87:     }
 88:   }
 89:   /* Process summary exclusions */
 90:   PetscOptionsGetString(NULL,NULL,"-log_exclude",logList,256,&opt);
 91:   if (opt) {
 92:     PetscStrstr(logList,"fn",&className);
 93:     if (className) {
 94:       PetscLogEventDeactivateClass(FN_CLASSID);
 95:     }
 96:   }
 97:   PetscRegisterFinalize(FNFinalizePackage);
 98:   return(0);
 99: }

103: /*@
104:    FNCreate - Creates an FN context.

106:    Collective on MPI_Comm

108:    Input Parameter:
109: .  comm - MPI communicator

111:    Output Parameter:
112: .  newfn - location to put the FN context

114:    Level: beginner

116: .seealso: FNDestroy(), FN
117: @*/
118: PetscErrorCode FNCreate(MPI_Comm comm,FN *newfn)
119: {
120:   FN             fn;

125:   *newfn = 0;
126:   FNInitializePackage();
127:   SlepcHeaderCreate(fn,FN_CLASSID,"FN","Math Function","FN",comm,FNDestroy,FNView);

129:   fn->alpha    = 1.0;
130:   fn->beta     = 1.0;

132:   fn->nw       = 0;
133:   fn->cw       = 0;
134:   fn->data     = NULL;

136:   *newfn = fn;
137:   return(0);
138: }

142: /*@C
143:    FNSetOptionsPrefix - Sets the prefix used for searching for all
144:    FN options in the database.

146:    Logically Collective on FN

148:    Input Parameters:
149: +  fn - the math function context
150: -  prefix - the prefix string to prepend to all FN option requests

152:    Notes:
153:    A hyphen (-) must NOT be given at the beginning of the prefix name.
154:    The first character of all runtime options is AUTOMATICALLY the
155:    hyphen.

157:    Level: advanced

159: .seealso: FNAppendOptionsPrefix()
160: @*/
161: PetscErrorCode FNSetOptionsPrefix(FN fn,const char *prefix)
162: {

167:   PetscObjectSetOptionsPrefix((PetscObject)fn,prefix);
168:   return(0);
169: }

173: /*@C
174:    FNAppendOptionsPrefix - Appends to the prefix used for searching for all
175:    FN options in the database.

177:    Logically Collective on FN

179:    Input Parameters:
180: +  fn - the math function context
181: -  prefix - the prefix string to prepend to all FN option requests

183:    Notes:
184:    A hyphen (-) must NOT be given at the beginning of the prefix name.
185:    The first character of all runtime options is AUTOMATICALLY the hyphen.

187:    Level: advanced

189: .seealso: FNSetOptionsPrefix()
190: @*/
191: PetscErrorCode FNAppendOptionsPrefix(FN fn,const char *prefix)
192: {

197:   PetscObjectAppendOptionsPrefix((PetscObject)fn,prefix);
198:   return(0);
199: }

203: /*@C
204:    FNGetOptionsPrefix - Gets the prefix used for searching for all
205:    FN options in the database.

207:    Not Collective

209:    Input Parameters:
210: .  fn - the math function context

212:    Output Parameters:
213: .  prefix - pointer to the prefix string used is returned

215:    Note:
216:    On the Fortran side, the user should pass in a string 'prefix' of
217:    sufficient length to hold the prefix.

219:    Level: advanced

221: .seealso: FNSetOptionsPrefix(), FNAppendOptionsPrefix()
222: @*/
223: PetscErrorCode FNGetOptionsPrefix(FN fn,const char *prefix[])
224: {

230:   PetscObjectGetOptionsPrefix((PetscObject)fn,prefix);
231:   return(0);
232: }

236: /*@C
237:    FNSetType - Selects the type for the FN object.

239:    Logically Collective on FN

241:    Input Parameter:
242: +  fn   - the math function context
243: -  type - a known type

245:    Notes:
246:    The default is FNRATIONAL, which includes polynomials as a particular
247:    case as well as simple functions such as f(x)=x and f(x)=constant.

249:    Level: intermediate

251: .seealso: FNGetType()
252: @*/
253: PetscErrorCode FNSetType(FN fn,FNType type)
254: {
255:   PetscErrorCode ierr,(*r)(FN);
256:   PetscBool      match;


262:   PetscObjectTypeCompare((PetscObject)fn,type,&match);
263:   if (match) return(0);

265:    PetscFunctionListFind(FNList,type,&r);
266:   if (!r) SETERRQ1(PetscObjectComm((PetscObject)fn),PETSC_ERR_ARG_UNKNOWN_TYPE,"Unable to find requested FN type %s",type);

268:   if (fn->ops->destroy) { (*fn->ops->destroy)(fn); }
269:   PetscMemzero(fn->ops,sizeof(struct _FNOps));

271:   PetscObjectChangeTypeName((PetscObject)fn,type);
272:   (*r)(fn);
273:   return(0);
274: }

278: /*@C
279:    FNGetType - Gets the FN type name (as a string) from the FN context.

281:    Not Collective

283:    Input Parameter:
284: .  fn - the math function context

286:    Output Parameter:
287: .  name - name of the math function

289:    Level: intermediate

291: .seealso: FNSetType()
292: @*/
293: PetscErrorCode FNGetType(FN fn,FNType *type)
294: {
298:   *type = ((PetscObject)fn)->type_name;
299:   return(0);
300: }

304: /*@
305:    FNSetScale - Sets the scaling parameters that define the matematical function.

307:    Logically Collective on FN

309:    Input Parameters:
310: +  fn    - the math function context
311: .  alpha - inner scaling (argument)
312: -  beta  - outer scaling (result)

314:    Notes:
315:    Given a function f(x) specified by the FN type, the scaling parameters can
316:    be used to realize the function beta*f(alpha*x). So when these values are given,
317:    the procedure for function evaluation will first multiply the argument by alpha,
318:    then evaluate the function itself, and finally scale the result by beta.
319:    Likewise, these values are also considered when evaluating the derivative.

321:    If you want to provide only one of the two scaling factors, set the other
322:    one to 1.0.

324:    Level: intermediate

326: .seealso: FNGetScale(), FNEvaluateFunction()
327: @*/
328: PetscErrorCode FNSetScale(FN fn,PetscScalar alpha,PetscScalar beta)
329: {
334:   fn->alpha = alpha;
335:   fn->beta  = beta;
336:   return(0);
337: }

341: /*@
342:    FNGetScale - Gets the scaling parameters that define the matematical function.

344:    Not Collective

346:    Input Parameter:
347: .  fn    - the math function context

349:    Output Parameters:
350: +  alpha - inner scaling (argument)
351: -  beta  - outer scaling (result)

353:    Level: intermediate

355: .seealso: FNSetScale()
356: @*/
357: PetscErrorCode FNGetScale(FN fn,PetscScalar *alpha,PetscScalar *beta)
358: {
361:   if (alpha) *alpha = fn->alpha;
362:   if (beta)  *beta  = fn->beta;
363:   return(0);
364: }

368: /*@
369:    FNEvaluateFunction - Computes the value of the function f(x) for a given x.

371:    Logically Collective on FN

373:    Input Parameters:
374: +  fn - the math function context
375: -  x  - the value where the function must be evaluated

377:    Output Parameter:
378: .  y  - the result of f(x)

380:    Note:
381:    Scaling factors are taken into account, so the actual function evaluation
382:    will return beta*f(alpha*x).

384:    Level: intermediate

386: .seealso: FNEvaluateDerivative(), FNEvaluateFunctionMat(), FNSetScale()
387: @*/
388: PetscErrorCode FNEvaluateFunction(FN fn,PetscScalar x,PetscScalar *y)
389: {
391:   PetscScalar    xf,yf;

398:   PetscLogEventBegin(FN_Evaluate,fn,0,0,0);
399:   xf = fn->alpha*x;
400:   (*fn->ops->evaluatefunction)(fn,xf,&yf);
401:   *y = fn->beta*yf;
402:   PetscLogEventEnd(FN_Evaluate,fn,0,0,0);
403:   return(0);
404: }

408: /*@
409:    FNEvaluateDerivative - Computes the value of the derivative f'(x) for a given x.

411:    Logically Collective on FN

413:    Input Parameters:
414: +  fn - the math function context
415: -  x  - the value where the derivative must be evaluated

417:    Output Parameter:
418: .  y  - the result of f'(x)

420:    Note:
421:    Scaling factors are taken into account, so the actual derivative evaluation will
422:    return alpha*beta*f'(alpha*x).

424:    Level: intermediate

426: .seealso: FNEvaluateFunction()
427: @*/
428: PetscErrorCode FNEvaluateDerivative(FN fn,PetscScalar x,PetscScalar *y)
429: {
431:   PetscScalar    xf,yf;

438:   PetscLogEventBegin(FN_Evaluate,fn,0,0,0);
439:   xf = fn->alpha*x;
440:   (*fn->ops->evaluatederivative)(fn,xf,&yf);
441:   *y = fn->alpha*fn->beta*yf;
442:   PetscLogEventEnd(FN_Evaluate,fn,0,0,0);
443:   return(0);
444: }

448: static PetscErrorCode FNEvaluateFunctionMat_Sym_Private(FN fn,PetscScalar *As,PetscScalar *Bs,PetscInt m,PetscBool firstonly)
449: {
450: #if defined(PETSC_MISSING_LAPACK_SYEV) || defined(SLEPC_MISSING_LAPACK_LACPY)
452:   SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"SYEV/LACPY - Lapack routines are unavailable");
453: #else
455:   PetscInt       i,j;
456:   PetscBLASInt   n,k,ld,lwork,info;
457:   PetscScalar    *Q,*W,*work,a,x,y,one=1.0,zero=0.0;
458:   PetscReal      *eig,dummy;
459: #if defined(PETSC_USE_COMPLEX)
460:   PetscReal      *rwork,rdummy;
461: #endif

464:   PetscBLASIntCast(m,&n);
465:   ld = n;
466:   k = firstonly? 1: n;

468:   /* workspace query and memory allocation */
469:   lwork = -1;
470: #if defined(PETSC_USE_COMPLEX)
471:   PetscStackCallBLAS("LAPACKsyev",LAPACKsyev_("V","L",&n,As,&ld,&dummy,&a,&lwork,&rdummy,&info));
472:   PetscBLASIntCast((PetscInt)PetscRealPart(a),&lwork);
473:   PetscMalloc5(m,&eig,m*m,&Q,m*k,&W,lwork,&work,PetscMax(1,3*m-2),&rwork);
474: #else
475:   PetscStackCallBLAS("LAPACKsyev",LAPACKsyev_("V","L",&n,As,&ld,&dummy,&a,&lwork,&info));
476:   PetscBLASIntCast((PetscInt)PetscRealPart(a),&lwork);
477:   PetscMalloc4(m,&eig,m*m,&Q,m*k,&W,lwork,&work);
478: #endif

480:   /* compute eigendecomposition */
481:   PetscStackCallBLAS("LAPACKlacpy",LAPACKlacpy_("L",&n,&n,As,&ld,Q,&ld));
482: #if defined(PETSC_USE_COMPLEX)
483:   PetscStackCallBLAS("LAPACKsyev",LAPACKsyev_("V","L",&n,Q,&ld,eig,work,&lwork,rwork,&info));
484: #else
485:   PetscStackCallBLAS("LAPACKsyev",LAPACKsyev_("V","L",&n,Q,&ld,eig,work,&lwork,&info));
486: #endif
487:   if (info) SETERRQ1(PetscObjectComm((PetscObject)fn),PETSC_ERR_LIB,"Error in Lapack xSYEV %i",info);

489:   /* W = f(Lambda)*Q' */
490:   for (i=0;i<n;i++) {
491:     x = eig[i];
492:     (*fn->ops->evaluatefunction)(fn,x,&y);  /* y = f(x) */
493:     for (j=0;j<k;j++) W[i+j*ld] = Q[j+i*ld]*y;
494:   }
495:   /* Bs = Q*W */
496:   PetscStackCallBLAS("BLASgemm",BLASgemm_("N","N",&n,&k,&n,&one,Q,&ld,W,&ld,&zero,Bs,&ld));
497: #if defined(PETSC_USE_COMPLEX)
498:   PetscFree5(eig,Q,W,work,rwork);
499: #else
500:   PetscFree4(eig,Q,W,work);
501: #endif
502:   return(0);
503: #endif
504: }

508: /*
509:    FNEvaluateFunctionMat_Sym_Default - given a symmetric matrix A,
510:    compute the matrix function as f(A)=Q*f(D)*Q' where the spectral
511:    decomposition of A is A=Q*D*Q'
512: */
513: static PetscErrorCode FNEvaluateFunctionMat_Sym_Default(FN fn,Mat A,Mat B)
514: {
516:   PetscInt       m;
517:   PetscScalar    *As,*Bs;

520:   MatDenseGetArray(A,&As);
521:   MatDenseGetArray(B,&Bs);
522:   MatGetSize(A,&m,NULL);
523:   FNEvaluateFunctionMat_Sym_Private(fn,As,Bs,m,PETSC_FALSE);
524:   MatDenseRestoreArray(A,&As);
525:   MatDenseRestoreArray(B,&Bs);
526:   return(0);
527: }

531: /*@
532:    FNEvaluateFunctionMat - Computes the value of the function f(A) for a given
533:    matrix A, where the result is also a matrix.

535:    Logically Collective on FN

537:    Input Parameters:
538: +  fn - the math function context
539: -  A  - matrix on which the function must be evaluated

541:    Output Parameter:
542: .  B  - (optional) matrix resulting from evaluating f(A)

544:    Notes:
545:    Matrix A must be a square sequential dense Mat, with all entries equal on
546:    all processes (otherwise each process will compute different results).
547:    If matrix B is provided, it must also be a square sequential dense Mat, and
548:    both matrices must have the same dimensions. If B is NULL (or B=A) then the
549:    function will perform an in-place computation, overwriting A with f(A).

551:    If A is known to be real symmetric or complex Hermitian then it is
552:    recommended to set the appropriate flag with MatSetOption(), so that
553:    a different algorithm that exploits symmetry is used.

555:    Scaling factors are taken into account, so the actual function evaluation
556:    will return beta*f(alpha*A).

558:    Level: advanced

560: .seealso: FNEvaluateFunction(), FNEvaluateFunctionMatVec()
561: @*/
562: PetscErrorCode FNEvaluateFunctionMat(FN fn,Mat A,Mat B)
563: {
565:   PetscBool      match,set,flg,symm=PETSC_FALSE,inplace=PETSC_FALSE;
566:   PetscInt       m,n,n1;
567:   Mat            M,F;

574:   if (B) {
577:   } else inplace = PETSC_TRUE;
578:   PetscObjectTypeCompare((PetscObject)A,MATSEQDENSE,&match);
579:   if (!match) SETERRQ(PetscObjectComm((PetscObject)fn),PETSC_ERR_SUP,"Mat A must be of type seqdense");
580:   MatGetSize(A,&m,&n);
581:   if (m!=n) SETERRQ2(PetscObjectComm((PetscObject)fn),PETSC_ERR_ARG_SIZ,"Mat A is not square (has %D rows, %D cols)",m,n);
582:   if (!inplace) {
583:     PetscObjectTypeCompare((PetscObject)B,MATSEQDENSE,&match);
584:     if (!match) SETERRQ(PetscObjectComm((PetscObject)fn),PETSC_ERR_SUP,"Mat B must be of type seqdense");
585:     n1 = n;
586:     MatGetSize(B,&m,&n);
587:     if (m!=n) SETERRQ2(PetscObjectComm((PetscObject)fn),PETSC_ERR_ARG_SIZ,"Mat B is not square (has %D rows, %D cols)",m,n);
588:     if (n1!=n) SETERRQ(PetscObjectComm((PetscObject)fn),PETSC_ERR_ARG_SIZ,"Matrices A and B must have the same dimension");
589:   }

591:   /* check symmetry of A */
592:   MatIsHermitianKnown(A,&set,&flg);
593:   symm = set? flg: PETSC_FALSE;

595:   /* scale argument */
596:   if (fn->alpha!=(PetscScalar)1.0) {
597:     FN_AllocateWorkMat(fn,A,&M);
598:     MatScale(M,fn->alpha);
599:   } else M = A;

601:   /* destination matrix */
602:   F = inplace? A: B;

604:   /* evaluate matrix function */
605:   PetscLogEventBegin(FN_Evaluate,fn,0,0,0);
606:   PetscFPTrapPush(PETSC_FP_TRAP_OFF);
607:   if (symm) {
608:     if (fn->ops->evaluatefunctionmatsym) {
609:       (*fn->ops->evaluatefunctionmatsym)(fn,M,F);
610:     } else {
611:       FNEvaluateFunctionMat_Sym_Default(fn,M,F);
612:     }
613:   } else {
614:     if (fn->ops->evaluatefunctionmat) {
615:       (*fn->ops->evaluatefunctionmat)(fn,M,F);
616:     } else SETERRQ1(PetscObjectComm((PetscObject)fn),PETSC_ERR_SUP,"Matrix function not implemented in FN type %s",((PetscObject)fn)->type_name);
617:   }
618:   PetscFPTrapPop();
619:   PetscLogEventEnd(FN_Evaluate,fn,0,0,0);

621:   if (fn->alpha!=(PetscScalar)1.0) {
622:     FN_FreeWorkMat(fn,&M);
623:   }

625:   /* scale result */
626:   MatScale(F,fn->beta);
627:   return(0);
628: }

632: /*
633:    FNEvaluateFunctionMatVec_Default - computes the full matrix f(A)
634:    and then copies the first column.
635: */
636: static PetscErrorCode FNEvaluateFunctionMatVec_Default(FN fn,Mat A,Vec v)
637: {
639:   Mat            F;

642:   FN_AllocateWorkMat(fn,A,&F);
643:   if (fn->ops->evaluatefunctionmat) {
644:     (*fn->ops->evaluatefunctionmat)(fn,A,F);
645:   } else SETERRQ1(PetscObjectComm((PetscObject)fn),PETSC_ERR_SUP,"Matrix function not implemented in FN type %s",((PetscObject)fn)->type_name);
646:   MatGetColumnVector(F,v,0);
647:   FN_FreeWorkMat(fn,&F);
648:   return(0);
649: }

653: /*
654:    FNEvaluateFunctionMatVec_Sym_Default - given a symmetric matrix A,
655:    compute the matrix function as f(A)=Q*f(D)*Q' where the spectral
656:    decomposition of A is A=Q*D*Q'. Only the first column is computed.
657: */
658: static PetscErrorCode FNEvaluateFunctionMatVec_Sym_Default(FN fn,Mat A,Vec v)
659: {
661:   PetscInt       m;
662:   PetscScalar    *As,*vs;

665:   MatDenseGetArray(A,&As);
666:   VecGetArray(v,&vs);
667:   MatGetSize(A,&m,NULL);
668:   FNEvaluateFunctionMat_Sym_Private(fn,As,vs,m,PETSC_TRUE);
669:   MatDenseRestoreArray(A,&As);
670:   VecRestoreArray(v,&vs);
671:   return(0);
672: }

676: /*@
677:    FNEvaluateFunctionMatVec - Computes the first column of the matrix f(A)
678:    for a given matrix A.

680:    Logically Collective on FN

682:    Input Parameters:
683: +  fn - the math function context
684: -  A  - matrix on which the function must be evaluated

686:    Output Parameter:
687: .  v  - vector to hold the first column of f(A)

689:    Notes:
690:    This operation is similar to FNEvaluateFunctionMat() but returns only
691:    the first column of f(A), hence saving computations in most cases.

693:    Level: advanced

695: .seealso: FNEvaluateFunction(), FNEvaluateFunctionMat()
696: @*/
697: PetscErrorCode FNEvaluateFunctionMatVec(FN fn,Mat A,Vec v)
698: {
700:   PetscBool      match,set,flg,symm=PETSC_FALSE;
701:   PetscInt       m,n;
702:   Mat            M;

711:   PetscObjectTypeCompare((PetscObject)A,MATSEQDENSE,&match);
712:   if (!match) SETERRQ(PetscObjectComm((PetscObject)fn),PETSC_ERR_SUP,"Mat A must be of type seqdense");
713:   MatGetSize(A,&m,&n);
714:   if (m!=n) SETERRQ2(PetscObjectComm((PetscObject)fn),PETSC_ERR_ARG_SIZ,"Mat A is not square (has %D rows, %D cols)",m,n);
715:   VecGetSize(v,&m);
716:   if (m!=n) SETERRQ(PetscObjectComm((PetscObject)fn),PETSC_ERR_ARG_SIZ,"Matrix A and vector v must have the same size");

718:   /* check symmetry of A */
719:   MatIsHermitianKnown(A,&set,&flg);
720:   symm = set? flg: PETSC_FALSE;

722:   /* scale argument */
723:   if (fn->alpha!=(PetscScalar)1.0) {
724:     FN_AllocateWorkMat(fn,A,&M);
725:     MatScale(M,fn->alpha);
726:   } else M = A;

728:   /* evaluate matrix function */
729:   PetscLogEventBegin(FN_Evaluate,fn,0,0,0);
730:   PetscFPTrapPush(PETSC_FP_TRAP_OFF);
731:   if (symm) {
732:     if (fn->ops->evaluatefunctionmatvecsym) {
733:       (*fn->ops->evaluatefunctionmatvecsym)(fn,M,v);
734:     } else {
735:       FNEvaluateFunctionMatVec_Sym_Default(fn,M,v);
736:     }
737:   } else {
738:     if (fn->ops->evaluatefunctionmatvec) {
739:       (*fn->ops->evaluatefunctionmatvec)(fn,M,v);
740:     } else {
741:       FNEvaluateFunctionMatVec_Default(fn,M,v);
742:     }
743:   }
744:   PetscFPTrapPop();
745:   PetscLogEventEnd(FN_Evaluate,fn,0,0,0);

747:   if (fn->alpha!=(PetscScalar)1.0) {
748:     FN_FreeWorkMat(fn,&M);
749:   }

751:   /* scale result */
752:   VecScale(v,fn->beta);
753:   return(0);
754: }

758: /*@
759:    FNSetFromOptions - Sets FN options from the options database.

761:    Collective on FN

763:    Input Parameters:
764: .  fn - the math function context

766:    Notes:
767:    To see all options, run your program with the -help option.

769:    Level: beginner
770: @*/
771: PetscErrorCode FNSetFromOptions(FN fn)
772: {
774:   char           type[256];
775:   PetscScalar    array[2];
776:   PetscInt       k;
777:   PetscBool      flg;

781:   FNRegisterAll();
782:   PetscObjectOptionsBegin((PetscObject)fn);
783:     PetscOptionsFList("-fn_type","Math function type","FNSetType",FNList,(char*)(((PetscObject)fn)->type_name?((PetscObject)fn)->type_name:FNRATIONAL),type,256,&flg);
784:     if (flg) {
785:       FNSetType(fn,type);
786:     }
787:     /*
788:       Set the type if it was never set.
789:     */
790:     if (!((PetscObject)fn)->type_name) {
791:       FNSetType(fn,FNRATIONAL);
792:     }

794:     k = 2;
795:     array[0] = 0.0; array[1] = 0.0;
796:     PetscOptionsScalarArray("-fn_scale","Scale factors (one or two scalar values separated with a comma without spaces)","FNSetScale",array,&k,&flg);
797:     if (flg) {
798:       if (k<2) array[1] = 1.0;
799:       FNSetScale(fn,array[0],array[1]);
800:     }

802:     if (fn->ops->setfromoptions) {
803:       (*fn->ops->setfromoptions)(PetscOptionsObject,fn);
804:     }
805:     PetscObjectProcessOptionsHandlers(PetscOptionsObject,(PetscObject)fn);
806:   PetscOptionsEnd();
807:   return(0);
808: }

812: /*@C
813:    FNView - Prints the FN data structure.

815:    Collective on FN

817:    Input Parameters:
818: +  fn - the math function context
819: -  viewer - optional visualization context

821:    Note:
822:    The available visualization contexts include
823: +     PETSC_VIEWER_STDOUT_SELF - standard output (default)
824: -     PETSC_VIEWER_STDOUT_WORLD - synchronized standard
825:          output where only the first processor opens
826:          the file.  All other processors send their
827:          data to the first processor to print.

829:    The user can open an alternative visualization context with
830:    PetscViewerASCIIOpen() - output to a specified file.

832:    Level: beginner
833: @*/
834: PetscErrorCode FNView(FN fn,PetscViewer viewer)
835: {
836:   PetscBool      isascii;

841:   if (!viewer) viewer = PETSC_VIEWER_STDOUT_(PetscObjectComm((PetscObject)fn));
844:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);
845:   if (isascii) {
846:     PetscObjectPrintClassNamePrefixType((PetscObject)fn,viewer);
847:     if (fn->ops->view) {
848:       PetscViewerASCIIPushTab(viewer);
849:       (*fn->ops->view)(fn,viewer);
850:       PetscViewerASCIIPopTab(viewer);
851:     }
852:   }
853:   return(0);
854: }

858: /*@
859:    FNDuplicate - Duplicates a math function, copying all parameters, possibly with a
860:    different communicator.

862:    Collective on FN

864:    Input Parameters:
865: +  fn   - the math function context
866: -  comm - MPI communicator

868:    Output Parameter:
869: .  newfn - location to put the new FN context

871:    Note:
872:    In order to use the same MPI communicator as in the original object,
873:    use PetscObjectComm((PetscObject)fn).

875:    Level: developer

877: .seealso: FNCreate()
878: @*/
879: PetscErrorCode FNDuplicate(FN fn,MPI_Comm comm,FN *newfn)
880: {
882:   FNType         type;
883:   PetscScalar    alpha,beta;

889:   FNCreate(comm,newfn);
890:   FNGetType(fn,&type);
891:   FNSetType(*newfn,type);
892:   FNGetScale(fn,&alpha,&beta);
893:   FNSetScale(*newfn,alpha,beta);
894:   if (fn->ops->duplicate) {
895:     (*fn->ops->duplicate)(fn,comm,newfn);
896:   }
897:   return(0);
898: }

902: /*@
903:    FNDestroy - Destroys FN context that was created with FNCreate().

905:    Collective on FN

907:    Input Parameter:
908: .  fn - the math function context

910:    Level: beginner

912: .seealso: FNCreate()
913: @*/
914: PetscErrorCode FNDestroy(FN *fn)
915: {
917:   PetscInt       i;

920:   if (!*fn) return(0);
922:   if (--((PetscObject)(*fn))->refct > 0) { *fn = 0; return(0); }
923:   if ((*fn)->ops->destroy) { (*(*fn)->ops->destroy)(*fn); }
924:   for (i=0;i<(*fn)->nw;i++) {
925:     MatDestroy(&(*fn)->W[i]);
926:   }
927:   PetscHeaderDestroy(fn);
928:   return(0);
929: }

933: /*@C
934:    FNRegister - Adds a mathematical function to the FN package.

936:    Not collective

938:    Input Parameters:
939: +  name - name of a new user-defined FN
940: -  function - routine to create context

942:    Notes:
943:    FNRegister() may be called multiple times to add several user-defined functions.

945:    Level: advanced

947: .seealso: FNRegisterAll()
948: @*/
949: PetscErrorCode FNRegister(const char *name,PetscErrorCode (*function)(FN))
950: {

954:   PetscFunctionListAdd(&FNList,name,function);
955:   return(0);
956: }