Actual source code: fnsqrt.c

slepc-3.7.0 2016-05-16
Report Typos and Errors
  1: /*
  2:    Square root function  sqrt(x)

  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>

 29: PetscErrorCode FNEvaluateFunction_Sqrt(FN fn,PetscScalar x,PetscScalar *y)
 30: {
 32:   *y = PetscSqrtScalar(x);
 33:   return(0);
 34: }

 38: PetscErrorCode FNEvaluateDerivative_Sqrt(FN fn,PetscScalar x,PetscScalar *y)
 39: {
 41:   if (x==0.0) SETERRQ(PETSC_COMM_SELF,1,"Derivative not defined in the requested value");
 42:   *y = 1.0/(2.0*PetscSqrtScalar(x));
 43:   return(0);
 44: }

 48: PetscErrorCode FNEvaluateFunctionMat_Sqrt(FN fn,Mat A,Mat B)
 49: {
 51:   PetscBLASInt   n;
 52:   PetscScalar    *T;
 53:   PetscInt       m;

 56:   if (A!=B) { MatCopy(A,B,SAME_NONZERO_PATTERN); }
 57:   MatDenseGetArray(B,&T);
 58:   MatGetSize(A,&m,NULL);
 59:   PetscBLASIntCast(m,&n);
 60:   SlepcSchurParlettSqrt(n,T,n,PETSC_FALSE);
 61:   MatDenseRestoreArray(B,&T);
 62:   return(0);
 63: }

 67: PetscErrorCode FNEvaluateFunctionMatVec_Sqrt(FN fn,Mat A,Vec v)
 68: {
 70:   PetscBLASInt   n;
 71:   PetscScalar    *T;
 72:   PetscInt       m;
 73:   Mat            B;

 76:   FN_AllocateWorkMat(fn,A,&B);
 77:   MatDenseGetArray(B,&T);
 78:   MatGetSize(A,&m,NULL);
 79:   PetscBLASIntCast(m,&n);
 80:   SlepcSchurParlettSqrt(n,T,n,PETSC_TRUE);
 81:   MatDenseRestoreArray(B,&T);
 82:   MatGetColumnVector(B,v,0);
 83:   FN_FreeWorkMat(fn,&B);
 84:   return(0);
 85: }

 89: PetscErrorCode FNView_Sqrt(FN fn,PetscViewer viewer)
 90: {
 92:   PetscBool      isascii;
 93:   char           str[50];

 96:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);
 97:   if (isascii) {
 98:     if (fn->beta==(PetscScalar)1.0) {
 99:       if (fn->alpha==(PetscScalar)1.0) {
100:         PetscViewerASCIIPrintf(viewer,"  Square root: sqrt(x)\n");
101:       } else {
102:         SlepcSNPrintfScalar(str,50,fn->alpha,PETSC_TRUE);
103:         PetscViewerASCIIPrintf(viewer,"  Square root: sqrt(%s*x)\n",str);
104:       }
105:     } else {
106:       SlepcSNPrintfScalar(str,50,fn->beta,PETSC_TRUE);
107:       if (fn->alpha==(PetscScalar)1.0) {
108:         PetscViewerASCIIPrintf(viewer,"  Square root: %s*sqrt(x)\n",str);
109:       } else {
110:         PetscViewerASCIIPrintf(viewer,"  Square root: %s",str);
111:         PetscViewerASCIIUseTabs(viewer,PETSC_FALSE);
112:         SlepcSNPrintfScalar(str,50,fn->alpha,PETSC_TRUE);
113:         PetscViewerASCIIPrintf(viewer,"*sqrt(%s*x)\n",str);
114:         PetscViewerASCIIUseTabs(viewer,PETSC_TRUE);
115:       }
116:     }
117:   }
118:   return(0);
119: }

123: PETSC_EXTERN PetscErrorCode FNCreate_Sqrt(FN fn)
124: {
126:   fn->ops->evaluatefunction       = FNEvaluateFunction_Sqrt;
127:   fn->ops->evaluatederivative     = FNEvaluateDerivative_Sqrt;
128:   fn->ops->evaluatefunctionmat    = FNEvaluateFunctionMat_Sqrt;
129:   fn->ops->evaluatefunctionmatvec = FNEvaluateFunctionMatVec_Sqrt;
130:   fn->ops->view                   = FNView_Sqrt;
131:   return(0);
132: }