Actual source code: fnbasic.c

slepc-3.8.2 2017-12-01
Report Typos and Errors
  1: /*
  2:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  3:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  4:    Copyright (c) 2002-2017, Universitat Politecnica de Valencia, Spain

  6:    This file is part of SLEPc.
  7:    SLEPc is distributed under a 2-clause BSD license (see LICENSE).
  8:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  9: */
 10: /*
 11:    Basic FN routines
 12: */

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

 17: PetscFunctionList FNList = 0;
 18: PetscBool         FNRegisterAllCalled = PETSC_FALSE;
 19: PetscClassId      FN_CLASSID = 0;
 20: PetscLogEvent     FN_Evaluate = 0;
 21: static PetscBool  FNPackageInitialized = PETSC_FALSE;

 23: const char *FNParallelTypes[] = {"REDUNDANT","SYNCHRONIZED","FNParallelType","FN_PARALLEL_",0};

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

 29:    Level: developer

 31: .seealso: SlepcFinalize()
 32: @*/
 33: PetscErrorCode FNFinalizePackage(void)
 34: {

 38:   PetscFunctionListDestroy(&FNList);
 39:   FNPackageInitialized = PETSC_FALSE;
 40:   FNRegisterAllCalled  = PETSC_FALSE;
 41:   return(0);
 42: }

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

 49:   Level: developer

 51: .seealso: SlepcInitialize()
 52: @*/
 53: PetscErrorCode FNInitializePackage(void)
 54: {
 55:   char             logList[256];
 56:   char             *className;
 57:   PetscBool        opt;
 58:   PetscErrorCode   ierr;

 61:   if (FNPackageInitialized) return(0);
 62:   FNPackageInitialized = PETSC_TRUE;
 63:   /* Register Classes */
 64:   PetscClassIdRegister("Math Function",&FN_CLASSID);
 65:   /* Register Constructors */
 66:   FNRegisterAll();
 67:   /* Register Events */
 68:   PetscLogEventRegister("FNEvaluate",FN_CLASSID,&FN_Evaluate);
 69:   /* Process info exclusions */
 70:   PetscOptionsGetString(NULL,NULL,"-info_exclude",logList,256,&opt);
 71:   if (opt) {
 72:     PetscStrstr(logList,"fn",&className);
 73:     if (className) {
 74:       PetscInfoDeactivateClass(FN_CLASSID);
 75:     }
 76:   }
 77:   /* Process summary exclusions */
 78:   PetscOptionsGetString(NULL,NULL,"-log_exclude",logList,256,&opt);
 79:   if (opt) {
 80:     PetscStrstr(logList,"fn",&className);
 81:     if (className) {
 82:       PetscLogEventDeactivateClass(FN_CLASSID);
 83:     }
 84:   }
 85:   PetscRegisterFinalize(FNFinalizePackage);
 86:   return(0);
 87: }

 89: /*@
 90:    FNCreate - Creates an FN context.

 92:    Collective on MPI_Comm

 94:    Input Parameter:
 95: .  comm - MPI communicator

 97:    Output Parameter:
 98: .  newfn - location to put the FN context

100:    Level: beginner

102: .seealso: FNDestroy(), FN
103: @*/
104: PetscErrorCode FNCreate(MPI_Comm comm,FN *newfn)
105: {
106:   FN             fn;

111:   *newfn = 0;
112:   FNInitializePackage();
113:   SlepcHeaderCreate(fn,FN_CLASSID,"FN","Math Function","FN",comm,FNDestroy,FNView);

115:   fn->alpha    = 1.0;
116:   fn->beta     = 1.0;
117:   fn->method   = 0;

119:   fn->nw       = 0;
120:   fn->cw       = 0;
121:   fn->data     = NULL;

123:   *newfn = fn;
124:   return(0);
125: }

127: /*@C
128:    FNSetOptionsPrefix - Sets the prefix used for searching for all
129:    FN options in the database.

131:    Logically Collective on FN

133:    Input Parameters:
134: +  fn - the math function context
135: -  prefix - the prefix string to prepend to all FN option requests

137:    Notes:
138:    A hyphen (-) must NOT be given at the beginning of the prefix name.
139:    The first character of all runtime options is AUTOMATICALLY the
140:    hyphen.

142:    Level: advanced

144: .seealso: FNAppendOptionsPrefix()
145: @*/
146: PetscErrorCode FNSetOptionsPrefix(FN fn,const char *prefix)
147: {

152:   PetscObjectSetOptionsPrefix((PetscObject)fn,prefix);
153:   return(0);
154: }

156: /*@C
157:    FNAppendOptionsPrefix - Appends to the prefix used for searching for all
158:    FN options in the database.

160:    Logically Collective on FN

162:    Input Parameters:
163: +  fn - the math function context
164: -  prefix - the prefix string to prepend to all FN option requests

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

170:    Level: advanced

172: .seealso: FNSetOptionsPrefix()
173: @*/
174: PetscErrorCode FNAppendOptionsPrefix(FN fn,const char *prefix)
175: {

180:   PetscObjectAppendOptionsPrefix((PetscObject)fn,prefix);
181:   return(0);
182: }

184: /*@C
185:    FNGetOptionsPrefix - Gets the prefix used for searching for all
186:    FN options in the database.

188:    Not Collective

190:    Input Parameters:
191: .  fn - the math function context

193:    Output Parameters:
194: .  prefix - pointer to the prefix string used is returned

196:    Note:
197:    On the Fortran side, the user should pass in a string 'prefix' of
198:    sufficient length to hold the prefix.

200:    Level: advanced

202: .seealso: FNSetOptionsPrefix(), FNAppendOptionsPrefix()
203: @*/
204: PetscErrorCode FNGetOptionsPrefix(FN fn,const char *prefix[])
205: {

211:   PetscObjectGetOptionsPrefix((PetscObject)fn,prefix);
212:   return(0);
213: }

215: /*@C
216:    FNSetType - Selects the type for the FN object.

218:    Logically Collective on FN

220:    Input Parameter:
221: +  fn   - the math function context
222: -  type - a known type

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

228:    Level: intermediate

230: .seealso: FNGetType()
231: @*/
232: PetscErrorCode FNSetType(FN fn,FNType type)
233: {
234:   PetscErrorCode ierr,(*r)(FN);
235:   PetscBool      match;


241:   PetscObjectTypeCompare((PetscObject)fn,type,&match);
242:   if (match) return(0);

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

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

250:   PetscObjectChangeTypeName((PetscObject)fn,type);
251:   (*r)(fn);
252:   return(0);
253: }

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

258:    Not Collective

260:    Input Parameter:
261: .  fn - the math function context

263:    Output Parameter:
264: .  name - name of the math function

266:    Level: intermediate

268: .seealso: FNSetType()
269: @*/
270: PetscErrorCode FNGetType(FN fn,FNType *type)
271: {
275:   *type = ((PetscObject)fn)->type_name;
276:   return(0);
277: }

279: /*@
280:    FNSetScale - Sets the scaling parameters that define the matematical function.

282:    Logically Collective on FN

284:    Input Parameters:
285: +  fn    - the math function context
286: .  alpha - inner scaling (argument)
287: -  beta  - outer scaling (result)

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

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

299:    Level: intermediate

301: .seealso: FNGetScale(), FNEvaluateFunction()
302: @*/
303: PetscErrorCode FNSetScale(FN fn,PetscScalar alpha,PetscScalar beta)
304: {
309:   fn->alpha = alpha;
310:   fn->beta  = beta;
311:   return(0);
312: }

314: /*@
315:    FNGetScale - Gets the scaling parameters that define the matematical function.

317:    Not Collective

319:    Input Parameter:
320: .  fn    - the math function context

322:    Output Parameters:
323: +  alpha - inner scaling (argument)
324: -  beta  - outer scaling (result)

326:    Level: intermediate

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

339: /*@
340:    FNSetMethod - Selects the method to be used to evaluate functions of matrices.

342:    Logically Collective on FN

344:    Input Parameter:
345: +  fn   - the math function context
346: -  meth - an index indentifying the method

348:    Options Database Key:
349: .  -fn_method <meth> - Sets the method

351:    Notes:
352:    In some FN types there are more than one algorithms available for computing
353:    matrix functions. In that case, this function allows choosing the wanted method.

355:    If meth is currently set to 0 (the default) and the input argument A of
356:    FNEvaluateFunctionMat() is a symmetric/Hermitian matrix, then the computation
357:    is done via the eigendecomposition of A, rather than with the general algorithm.

359:    Level: intermediate

361: .seealso: FNGetMethod(), FNEvaluateFunctionMat()
362: @*/
363: PetscErrorCode FNSetMethod(FN fn,PetscInt meth)
364: {
368:   if (meth<0) SETERRQ(PetscObjectComm((PetscObject)fn),PETSC_ERR_ARG_OUTOFRANGE,"The method must be a non-negative integer");
369:   if (meth>FN_MAX_SOLVE) SETERRQ(PetscObjectComm((PetscObject)fn),PETSC_ERR_ARG_OUTOFRANGE,"Too large value for the method");
370:   fn->method = meth;
371:   return(0);
372: }

374: /*@
375:    FNGetMethod - Gets the method currently used in the FN.

377:    Not Collective

379:    Input Parameter:
380: .  fn - the math function context

382:    Output Parameter:
383: .  meth - identifier of the method

385:    Level: intermediate

387: .seealso: FNSetMethod()
388: @*/
389: PetscErrorCode FNGetMethod(FN fn,PetscInt *meth)
390: {
394:   *meth = fn->method;
395:   return(0);
396: }

398: /*@
399:    FNSetParallel - Selects the mode of operation in parallel runs.

401:    Logically Collective on FN

403:    Input Parameter:
404: +  fn    - the math function context
405: -  pmode - the parallel mode

407:    Options Database Key:
408: .  -fn_parallel <mode> - Sets the parallel mode, either 'redundant' or 'synchronized'

410:    Notes:
411:    This is relevant only when the function is evaluated on a matrix, with
412:    either FNEvaluateFunctionMat() or FNEvaluateFunctionMatVec().

414:    In the 'redundant' parallel mode, all processes will make the computation
415:    redundantly, starting from the same data, and producing the same result.
416:    This result may be slightly different in the different processes if using a
417:    multithreaded BLAS library, which may cause issues in ill-conditioned problems.

419:    In the 'synchronized' parallel mode, only the first MPI process performs the
420:    computation and then the computed matrix is broadcast to the other
421:    processes in the communicator. This communication is done automatically at
422:    the end of FNEvaluateFunctionMat() or FNEvaluateFunctionMatVec().

424:    Level: advanced

426: .seealso: FNEvaluateFunctionMat() or FNEvaluateFunctionMatVec(), FNGetParallel()
427: @*/
428: PetscErrorCode FNSetParallel(FN fn,FNParallelType pmode)
429: {
433:   fn->pmode = pmode;
434:   return(0);
435: }

437: /*@
438:    FNGetParallel - Gets the mode of operation in parallel runs.

440:    Not Collective

442:    Input Parameter:
443: .  fn - the math function context

445:    Output Parameter:
446: .  pmode - the parallel mode

448:    Level: advanced

450: .seealso: FNSetParallel()
451: @*/
452: PetscErrorCode FNGetParallel(FN fn,FNParallelType *pmode)
453: {
457:   *pmode = fn->pmode;
458:   return(0);
459: }

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

464:    Not collective

466:    Input Parameters:
467: +  fn - the math function context
468: -  x  - the value where the function must be evaluated

470:    Output Parameter:
471: .  y  - the result of f(x)

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

477:    Level: intermediate

479: .seealso: FNEvaluateDerivative(), FNEvaluateFunctionMat(), FNSetScale()
480: @*/
481: PetscErrorCode FNEvaluateFunction(FN fn,PetscScalar x,PetscScalar *y)
482: {
484:   PetscScalar    xf,yf;

490:   PetscLogEventBegin(FN_Evaluate,fn,0,0,0);
491:   xf = fn->alpha*x;
492:   (*fn->ops->evaluatefunction)(fn,xf,&yf);
493:   *y = fn->beta*yf;
494:   PetscLogEventEnd(FN_Evaluate,fn,0,0,0);
495:   return(0);
496: }

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

501:    Not Collective

503:    Input Parameters:
504: +  fn - the math function context
505: -  x  - the value where the derivative must be evaluated

507:    Output Parameter:
508: .  y  - the result of f'(x)

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

514:    Level: intermediate

516: .seealso: FNEvaluateFunction()
517: @*/
518: PetscErrorCode FNEvaluateDerivative(FN fn,PetscScalar x,PetscScalar *y)
519: {
521:   PetscScalar    xf,yf;

527:   PetscLogEventBegin(FN_Evaluate,fn,0,0,0);
528:   xf = fn->alpha*x;
529:   (*fn->ops->evaluatederivative)(fn,xf,&yf);
530:   *y = fn->alpha*fn->beta*yf;
531:   PetscLogEventEnd(FN_Evaluate,fn,0,0,0);
532:   return(0);
533: }

535: static PetscErrorCode FNEvaluateFunctionMat_Sym_Private(FN fn,PetscScalar *As,PetscScalar *Bs,PetscInt m,PetscBool firstonly)
536: {
537: #if defined(PETSC_MISSING_LAPACK_SYEV) || defined(SLEPC_MISSING_LAPACK_LACPY)
539:   SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"SYEV/LACPY - Lapack routines are unavailable");
540: #else
542:   PetscInt       i,j;
543:   PetscBLASInt   n,k,ld,lwork,info;
544:   PetscScalar    *Q,*W,*work,a,x,y,one=1.0,zero=0.0;
545:   PetscReal      *eig,dummy;
546: #if defined(PETSC_USE_COMPLEX)
547:   PetscReal      *rwork,rdummy;
548: #endif

551:   PetscBLASIntCast(m,&n);
552:   ld = n;
553:   k = firstonly? 1: n;

555:   /* workspace query and memory allocation */
556:   lwork = -1;
557: #if defined(PETSC_USE_COMPLEX)
558:   PetscStackCallBLAS("LAPACKsyev",LAPACKsyev_("V","L",&n,As,&ld,&dummy,&a,&lwork,&rdummy,&info));
559:   PetscBLASIntCast((PetscInt)PetscRealPart(a),&lwork);
560:   PetscMalloc5(m,&eig,m*m,&Q,m*k,&W,lwork,&work,PetscMax(1,3*m-2),&rwork);
561: #else
562:   PetscStackCallBLAS("LAPACKsyev",LAPACKsyev_("V","L",&n,As,&ld,&dummy,&a,&lwork,&info));
563:   PetscBLASIntCast((PetscInt)PetscRealPart(a),&lwork);
564:   PetscMalloc4(m,&eig,m*m,&Q,m*k,&W,lwork,&work);
565: #endif

567:   /* compute eigendecomposition */
568:   PetscStackCallBLAS("LAPACKlacpy",LAPACKlacpy_("L",&n,&n,As,&ld,Q,&ld));
569: #if defined(PETSC_USE_COMPLEX)
570:   PetscStackCallBLAS("LAPACKsyev",LAPACKsyev_("V","L",&n,Q,&ld,eig,work,&lwork,rwork,&info));
571: #else
572:   PetscStackCallBLAS("LAPACKsyev",LAPACKsyev_("V","L",&n,Q,&ld,eig,work,&lwork,&info));
573: #endif
574:   SlepcCheckLapackInfo("syev",info);

576:   /* W = f(Lambda)*Q' */
577:   for (i=0;i<n;i++) {
578:     x = fn->alpha*eig[i];
579:     (*fn->ops->evaluatefunction)(fn,x,&y);  /* y = f(x) */
580:     for (j=0;j<k;j++) W[i+j*ld] = PetscConj(Q[j+i*ld])*fn->beta*y;
581:   }
582:   /* Bs = Q*W */
583:   PetscStackCallBLAS("BLASgemm",BLASgemm_("N","N",&n,&k,&n,&one,Q,&ld,W,&ld,&zero,Bs,&ld));
584: #if defined(PETSC_USE_COMPLEX)
585:   PetscFree5(eig,Q,W,work,rwork);
586: #else
587:   PetscFree4(eig,Q,W,work);
588: #endif
589:   PetscLogFlops(9.0*n*n*n+2.0*n*n*n);
590:   return(0);
591: #endif
592: }

594: /*
595:    FNEvaluateFunctionMat_Sym_Default - given a symmetric matrix A,
596:    compute the matrix function as f(A)=Q*f(D)*Q' where the spectral
597:    decomposition of A is A=Q*D*Q'
598: */
599: static PetscErrorCode FNEvaluateFunctionMat_Sym_Default(FN fn,Mat A,Mat B)
600: {
602:   PetscInt       m;
603:   PetscScalar    *As,*Bs;

606:   MatDenseGetArray(A,&As);
607:   MatDenseGetArray(B,&Bs);
608:   MatGetSize(A,&m,NULL);
609:   FNEvaluateFunctionMat_Sym_Private(fn,As,Bs,m,PETSC_FALSE);
610:   MatDenseRestoreArray(A,&As);
611:   MatDenseRestoreArray(B,&Bs);
612:   return(0);
613: }

615: PetscErrorCode FNEvaluateFunctionMat_Private(FN fn,Mat A,Mat B,PetscBool sync)
616: {
618:   PetscBool      set,flg,symm=PETSC_FALSE;
619:   PetscInt       m,n;
620:   PetscMPIInt    size,rank;
621:   PetscScalar    *pF;
622:   Mat            M,F;

625:   /* destination matrix */
626:   F = B?B:A;

628:   /* check symmetry of A */
629:   MatIsHermitianKnown(A,&set,&flg);
630:   symm = set? flg: PETSC_FALSE;

632:   MPI_Comm_size(PetscObjectComm((PetscObject)fn),&size);
633:   MPI_Comm_rank(PetscObjectComm((PetscObject)fn),&rank);
634:   if (size==1 || fn->pmode==FN_PARALLEL_REDUNDANT || (fn->pmode==FN_PARALLEL_SYNCHRONIZED && !rank)) {
635:     
636:     PetscFPTrapPush(PETSC_FP_TRAP_OFF);
637:     if (symm && !fn->method) {  /* prefer diagonalization */
638:       PetscInfo(fn,"Computing matrix function via diagonalization\n");
639:       FNEvaluateFunctionMat_Sym_Default(fn,A,F);
640:     } else {
641:       /* scale argument */
642:       if (fn->alpha!=(PetscScalar)1.0) {
643:         FN_AllocateWorkMat(fn,A,&M);
644:         MatScale(M,fn->alpha);
645:       } else M = A;
646:       if (fn->ops->evaluatefunctionmat[fn->method]) {
647:         (*fn->ops->evaluatefunctionmat[fn->method])(fn,M,F);
648:       } else SETERRQ(PetscObjectComm((PetscObject)fn),PETSC_ERR_ARG_OUTOFRANGE,"The specified method number does not exist for this FN");
649:       if (fn->alpha!=(PetscScalar)1.0) {
650:         FN_FreeWorkMat(fn,&M);
651:       }
652:       /* scale result */
653:       MatScale(F,fn->beta);
654:     }
655:     PetscFPTrapPop();
656:   }
657:   if (size>1 && fn->pmode==FN_PARALLEL_SYNCHRONIZED && sync) {  /* synchronize */
658:     MatGetSize(A,&m,&n);
659:     MatDenseGetArray(F,&pF);
660:     MPI_Bcast(pF,n*n,MPIU_SCALAR,0,PetscObjectComm((PetscObject)fn));
661:     MatDenseRestoreArray(F,&pF);
662:   }
663:   return(0);
664: }

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

670:    Logically Collective on FN

672:    Input Parameters:
673: +  fn - the math function context
674: -  A  - matrix on which the function must be evaluated

676:    Output Parameter:
677: .  B  - (optional) matrix resulting from evaluating f(A)

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

686:    If A is known to be real symmetric or complex Hermitian then it is
687:    recommended to set the appropriate flag with MatSetOption(), because
688:    symmetry can sometimes be exploited by the algorithm.

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

693:    Level: advanced

695: .seealso: FNEvaluateFunction(), FNEvaluateFunctionMatVec(), FNSetMethod()
696: @*/
697: PetscErrorCode FNEvaluateFunctionMat(FN fn,Mat A,Mat B)
698: {
700:   PetscBool      match,inplace=PETSC_FALSE;
701:   PetscInt       m,n,n1;

708:   if (B) {
711:   } else inplace = PETSC_TRUE;
712:   PetscObjectTypeCompare((PetscObject)A,MATSEQDENSE,&match);
713:   if (!match) SETERRQ(PetscObjectComm((PetscObject)fn),PETSC_ERR_SUP,"Mat A must be of type seqdense");
714:   MatGetSize(A,&m,&n);
715:   if (m!=n) SETERRQ2(PetscObjectComm((PetscObject)fn),PETSC_ERR_ARG_SIZ,"Mat A is not square (has %D rows, %D cols)",m,n);
716:   if (!inplace) {
717:     PetscObjectTypeCompare((PetscObject)B,MATSEQDENSE,&match);
718:     if (!match) SETERRQ(PetscObjectComm((PetscObject)fn),PETSC_ERR_SUP,"Mat B must be of type seqdense");
719:     n1 = n;
720:     MatGetSize(B,&m,&n);
721:     if (m!=n) SETERRQ2(PetscObjectComm((PetscObject)fn),PETSC_ERR_ARG_SIZ,"Mat B is not square (has %D rows, %D cols)",m,n);
722:     if (n1!=n) SETERRQ(PetscObjectComm((PetscObject)fn),PETSC_ERR_ARG_SIZ,"Matrices A and B must have the same dimension");
723:   }

725:   /* evaluate matrix function */
726:   PetscLogEventBegin(FN_Evaluate,fn,0,0,0);
727:   FNEvaluateFunctionMat_Private(fn,A,B,PETSC_TRUE);
728:   PetscLogEventEnd(FN_Evaluate,fn,0,0,0);
729:   return(0);
730: }

732: /*
733:    FNEvaluateFunctionMatVec_Default - computes the full matrix f(A)
734:    and then copies the first column.
735: */
736: static PetscErrorCode FNEvaluateFunctionMatVec_Default(FN fn,Mat A,Vec v)
737: {
739:   Mat            F;

742:   FN_AllocateWorkMat(fn,A,&F);
743:   if (fn->ops->evaluatefunctionmat[fn->method]) {
744:     (*fn->ops->evaluatefunctionmat[fn->method])(fn,A,F);
745:   } else SETERRQ(PetscObjectComm((PetscObject)fn),PETSC_ERR_ARG_OUTOFRANGE,"The specified method number does not exist for this FN");
746:   MatGetColumnVector(F,v,0);
747:   FN_FreeWorkMat(fn,&F);
748:   return(0);
749: }

751: /*
752:    FNEvaluateFunctionMatVec_Sym_Default - given a symmetric matrix A,
753:    compute the matrix function as f(A)=Q*f(D)*Q' where the spectral
754:    decomposition of A is A=Q*D*Q'. Only the first column is computed.
755: */
756: static PetscErrorCode FNEvaluateFunctionMatVec_Sym_Default(FN fn,Mat A,Vec v)
757: {
759:   PetscInt       m;
760:   PetscScalar    *As,*vs;

763:   MatDenseGetArray(A,&As);
764:   VecGetArray(v,&vs);
765:   MatGetSize(A,&m,NULL);
766:   FNEvaluateFunctionMat_Sym_Private(fn,As,vs,m,PETSC_TRUE);
767:   MatDenseRestoreArray(A,&As);
768:   VecRestoreArray(v,&vs);
769:   return(0);
770: }

772: PetscErrorCode FNEvaluateFunctionMatVec_Private(FN fn,Mat A,Vec v,PetscBool sync)
773: {
775:   PetscBool      set,flg,symm=PETSC_FALSE;
776:   PetscInt       m,n;
777:   Mat            M;
778:   PetscMPIInt    size,rank;
779:   PetscScalar    *pv;

782:   /* check symmetry of A */
783:   MatIsHermitianKnown(A,&set,&flg);
784:   symm = set? flg: PETSC_FALSE;

786:   /* evaluate matrix function */
787:   MPI_Comm_size(PetscObjectComm((PetscObject)fn),&size);
788:   MPI_Comm_rank(PetscObjectComm((PetscObject)fn),&rank);
789:   if (size==1 || fn->pmode==FN_PARALLEL_REDUNDANT || (fn->pmode==FN_PARALLEL_SYNCHRONIZED && !rank)) {
790:     PetscFPTrapPush(PETSC_FP_TRAP_OFF);
791:     if (symm && !fn->method) {  /* prefer diagonalization */
792:       PetscInfo(fn,"Computing matrix function via diagonalization\n");
793:       FNEvaluateFunctionMatVec_Sym_Default(fn,A,v);
794:     } else {
795:       /* scale argument */
796:       if (fn->alpha!=(PetscScalar)1.0) {
797:         FN_AllocateWorkMat(fn,A,&M);
798:         MatScale(M,fn->alpha);
799:       } else M = A;
800:       if (fn->ops->evaluatefunctionmatvec[fn->method]) {
801:         (*fn->ops->evaluatefunctionmatvec[fn->method])(fn,M,v);
802:       } else {
803:         FNEvaluateFunctionMatVec_Default(fn,M,v);
804:       }
805:       if (fn->alpha!=(PetscScalar)1.0) {
806:         FN_FreeWorkMat(fn,&M);
807:       }
808:       /* scale result */
809:       VecScale(v,fn->beta);
810:     }
811:     PetscFPTrapPop();
812:   }

814:   /* synchronize */
815:   if (size>1 && fn->pmode==FN_PARALLEL_SYNCHRONIZED && sync) { 
816:     MatGetSize(A,&m,&n);
817:     VecGetArray(v,&pv);
818:     MPI_Bcast(pv,n,MPIU_SCALAR,0,PetscObjectComm((PetscObject)fn));
819:     VecRestoreArray(v,&pv);
820:   }
821:   return(0);
822: }

824: /*@
825:    FNEvaluateFunctionMatVec - Computes the first column of the matrix f(A)
826:    for a given matrix A.

828:    Logically Collective on FN

830:    Input Parameters:
831: +  fn - the math function context
832: -  A  - matrix on which the function must be evaluated

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

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

841:    Level: advanced

843: .seealso: FNEvaluateFunction(), FNEvaluateFunctionMat(), FNSetMethod()
844: @*/
845: PetscErrorCode FNEvaluateFunctionMatVec(FN fn,Mat A,Vec v)
846: {
848:   PetscBool      match;
849:   PetscInt       m,n;

858:   PetscObjectTypeCompare((PetscObject)A,MATSEQDENSE,&match);
859:   if (!match) SETERRQ(PetscObjectComm((PetscObject)fn),PETSC_ERR_SUP,"Mat A must be of type seqdense");
860:   MatGetSize(A,&m,&n);
861:   if (m!=n) SETERRQ2(PetscObjectComm((PetscObject)fn),PETSC_ERR_ARG_SIZ,"Mat A is not square (has %D rows, %D cols)",m,n);
862:   VecGetSize(v,&m);
863:   if (m!=n) SETERRQ(PetscObjectComm((PetscObject)fn),PETSC_ERR_ARG_SIZ,"Matrix A and vector v must have the same size");
864:   PetscLogEventBegin(FN_Evaluate,fn,0,0,0);
865:   FNEvaluateFunctionMatVec_Private(fn,A,v,PETSC_TRUE);
866:   PetscLogEventEnd(FN_Evaluate,fn,0,0,0);
867:   return(0);
868: }

870: /*@
871:    FNSetFromOptions - Sets FN options from the options database.

873:    Collective on FN

875:    Input Parameters:
876: .  fn - the math function context

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

881:    Level: beginner
882: @*/
883: PetscErrorCode FNSetFromOptions(FN fn)
884: {
886:   char           type[256];
887:   PetscScalar    array[2];
888:   PetscInt       k,meth;
889:   PetscBool      flg;
890:   FNParallelType pmode;

894:   FNRegisterAll();
895:   PetscObjectOptionsBegin((PetscObject)fn);
896:     PetscOptionsFList("-fn_type","Math function type","FNSetType",FNList,(char*)(((PetscObject)fn)->type_name?((PetscObject)fn)->type_name:FNRATIONAL),type,256,&flg);
897:     if (flg) {
898:       FNSetType(fn,type);
899:     } else if (!((PetscObject)fn)->type_name) {
900:       FNSetType(fn,FNRATIONAL);
901:     }

903:     k = 2;
904:     array[0] = 0.0; array[1] = 0.0;
905:     PetscOptionsScalarArray("-fn_scale","Scale factors (one or two scalar values separated with a comma without spaces)","FNSetScale",array,&k,&flg);
906:     if (flg) {
907:       if (k<2) array[1] = 1.0;
908:       FNSetScale(fn,array[0],array[1]);
909:     }

911:     PetscOptionsInt("-fn_method","Method to be used for computing matrix functions","FNSetMethod",fn->method,&meth,&flg);
912:     if (flg) { FNSetMethod(fn,meth); }

914:     PetscOptionsEnum("-fn_parallel","Operation mode in parallel runs","FNSetParallel",FNParallelTypes,(PetscEnum)fn->pmode,(PetscEnum*)&pmode,&flg);
915:     if (flg) { FNSetParallel(fn,pmode); }

917:     if (fn->ops->setfromoptions) {
918:       (*fn->ops->setfromoptions)(PetscOptionsObject,fn);
919:     }
920:     PetscObjectProcessOptionsHandlers(PetscOptionsObject,(PetscObject)fn);
921:   PetscOptionsEnd();
922:   return(0);
923: }

925: /*@C
926:    FNView - Prints the FN data structure.

928:    Collective on FN

930:    Input Parameters:
931: +  fn - the math function context
932: -  viewer - optional visualization context

934:    Note:
935:    The available visualization contexts include
936: +     PETSC_VIEWER_STDOUT_SELF - standard output (default)
937: -     PETSC_VIEWER_STDOUT_WORLD - synchronized standard
938:          output where only the first processor opens
939:          the file.  All other processors send their
940:          data to the first processor to print.

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

945:    Level: beginner
946: @*/
947: PetscErrorCode FNView(FN fn,PetscViewer viewer)
948: {
949:   PetscBool      isascii;
951:   PetscMPIInt    size;

955:   if (!viewer) viewer = PETSC_VIEWER_STDOUT_(PetscObjectComm((PetscObject)fn));
958:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);
959:   if (isascii) {
960:     PetscObjectPrintClassNamePrefixType((PetscObject)fn,viewer);
961:     MPI_Comm_size(PetscObjectComm((PetscObject)fn),&size);
962:     if (size>1) {
963:       PetscViewerASCIIPrintf(viewer,"  parallel operation mode: %s\n",FNParallelTypes[fn->pmode]);
964:     }
965:     if (fn->ops->view) {
966:       PetscViewerASCIIPushTab(viewer);
967:       (*fn->ops->view)(fn,viewer);
968:       PetscViewerASCIIPopTab(viewer);
969:     }
970:   }
971:   return(0);
972: }

974: /*@
975:    FNDuplicate - Duplicates a math function, copying all parameters, possibly with a
976:    different communicator.

978:    Collective on FN

980:    Input Parameters:
981: +  fn   - the math function context
982: -  comm - MPI communicator

984:    Output Parameter:
985: .  newfn - location to put the new FN context

987:    Note:
988:    In order to use the same MPI communicator as in the original object,
989:    use PetscObjectComm((PetscObject)fn).

991:    Level: developer

993: .seealso: FNCreate()
994: @*/
995: PetscErrorCode FNDuplicate(FN fn,MPI_Comm comm,FN *newfn)
996: {
998:   FNType         type;
999:   PetscScalar    alpha,beta;
1000:   PetscInt       meth;
1001:   FNParallelType ptype;

1007:   FNCreate(comm,newfn);
1008:   FNGetType(fn,&type);
1009:   FNSetType(*newfn,type);
1010:   FNGetScale(fn,&alpha,&beta);
1011:   FNSetScale(*newfn,alpha,beta);
1012:   FNGetMethod(fn,&meth);
1013:   FNSetMethod(*newfn,meth);
1014:   FNGetParallel(fn,&ptype);
1015:   FNSetParallel(*newfn,ptype);
1016:   if (fn->ops->duplicate) {
1017:     (*fn->ops->duplicate)(fn,comm,newfn);
1018:   }
1019:   return(0);
1020: }

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

1025:    Collective on FN

1027:    Input Parameter:
1028: .  fn - the math function context

1030:    Level: beginner

1032: .seealso: FNCreate()
1033: @*/
1034: PetscErrorCode FNDestroy(FN *fn)
1035: {
1037:   PetscInt       i;

1040:   if (!*fn) return(0);
1042:   if (--((PetscObject)(*fn))->refct > 0) { *fn = 0; return(0); }
1043:   if ((*fn)->ops->destroy) { (*(*fn)->ops->destroy)(*fn); }
1044:   for (i=0;i<(*fn)->nw;i++) {
1045:     MatDestroy(&(*fn)->W[i]);
1046:   }
1047:   PetscHeaderDestroy(fn);
1048:   return(0);
1049: }

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

1054:    Not collective

1056:    Input Parameters:
1057: +  name - name of a new user-defined FN
1058: -  function - routine to create context

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

1063:    Level: advanced

1065: .seealso: FNRegisterAll()
1066: @*/
1067: PetscErrorCode FNRegister(const char *name,PetscErrorCode (*function)(FN))
1068: {

1072:   PetscFunctionListAdd(&FNList,name,function);
1073:   return(0);
1074: }