Actual source code: fninvsqrt.c

slepc-3.7.3 2016-09-29
Report Typos and Errors
  1: /*
  2:    Inverse square root function  x^(-1/2)

  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_Invsqrt(FN fn,PetscScalar x,PetscScalar *y)
 30: {
 32:   if (x==0.0) SETERRQ(PETSC_COMM_SELF,1,"Function not defined in the requested value");
 33:   *y = 1.0/PetscSqrtScalar(x);
 34:   return(0);
 35: }

 39: PetscErrorCode FNEvaluateDerivative_Invsqrt(FN fn,PetscScalar x,PetscScalar *y)
 40: {
 42:   if (x==0.0) SETERRQ(PETSC_COMM_SELF,1,"Derivative not defined in the requested value");
 43:   *y = -1.0/(2.0*PetscPowScalarReal(x,1.5));
 44:   return(0);
 45: }

 49: PetscErrorCode FNEvaluateFunctionMat_Invsqrt(FN fn,Mat A,Mat B)
 50: {
 52:   PetscBLASInt   n,ld,*ipiv,info;
 53:   PetscScalar    *Ba,*Wa;
 54:   PetscInt       m;
 55:   Mat            W;

 58:   FN_AllocateWorkMat(fn,A,&W);
 59:   if (A!=B) { MatCopy(A,B,SAME_NONZERO_PATTERN); }
 60:   MatDenseGetArray(B,&Ba);
 61:   MatDenseGetArray(W,&Wa);
 62:   /* compute B = sqrtm(A) */
 63:   MatGetSize(A,&m,NULL);
 64:   PetscBLASIntCast(m,&n);
 65:   ld = n;
 66:   SlepcSchurParlettSqrt(n,Ba,n,PETSC_FALSE);
 67:   /* compute B = A\B */
 68:   PetscMalloc1(ld,&ipiv);
 69:   PetscStackCallBLAS("LAPACKgesv",LAPACKgesv_(&n,&n,Wa,&ld,ipiv,Ba,&ld,&info));
 70:   if (info) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in Lapack xGESV %d",info);
 71:   PetscFree(ipiv);
 72:   MatDenseRestoreArray(W,&Wa);
 73:   MatDenseRestoreArray(B,&Ba);
 74:   FN_FreeWorkMat(fn,&W);
 75:   return(0);
 76: }

 80: PetscErrorCode FNEvaluateFunctionMatVec_Invsqrt(FN fn,Mat A,Vec v)
 81: {
 83:   PetscBLASInt   n,ld,*ipiv,info,one=1;
 84:   PetscScalar    *Ba,*Wa;
 85:   PetscInt       m;
 86:   Mat            B,W;

 89:   FN_AllocateWorkMat(fn,A,&B);
 90:   FN_AllocateWorkMat(fn,A,&W);
 91:   MatDenseGetArray(B,&Ba);
 92:   MatDenseGetArray(W,&Wa);
 93:   /* compute B_1 = sqrtm(A)*e_1 */
 94:   MatGetSize(A,&m,NULL);
 95:   PetscBLASIntCast(m,&n);
 96:   ld = n;
 97:   SlepcSchurParlettSqrt(n,Ba,n,PETSC_TRUE);
 98:   /* compute B_1 = A\B_1 */
 99:   PetscMalloc1(ld,&ipiv);
100:   PetscStackCallBLAS("LAPACKgesv",LAPACKgesv_(&n,&one,Wa,&ld,ipiv,Ba,&ld,&info));
101:   if (info) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in Lapack xGESV %d",info);
102:   PetscFree(ipiv);
103:   MatDenseRestoreArray(W,&Wa);
104:   MatDenseRestoreArray(B,&Ba);
105:   MatGetColumnVector(B,v,0);
106:   FN_FreeWorkMat(fn,&W);
107:   FN_FreeWorkMat(fn,&B);
108:   return(0);
109: }

113: PetscErrorCode FNView_Invsqrt(FN fn,PetscViewer viewer)
114: {
116:   PetscBool      isascii;
117:   char           str[50];

120:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);
121:   if (isascii) {
122:     if (fn->beta==(PetscScalar)1.0) {
123:       if (fn->alpha==(PetscScalar)1.0) {
124:         PetscViewerASCIIPrintf(viewer,"  Inverse square root: x^(-1/2)\n");
125:       } else {
126:         SlepcSNPrintfScalar(str,50,fn->alpha,PETSC_TRUE);
127:         PetscViewerASCIIPrintf(viewer,"  Inverse square root: (%s*x)^(-1/2)\n",str);
128:       }
129:     } else {
130:       SlepcSNPrintfScalar(str,50,fn->beta,PETSC_TRUE);
131:       if (fn->alpha==(PetscScalar)1.0) {
132:         PetscViewerASCIIPrintf(viewer,"  Inverse square root: %s*x^(-1/2)\n",str);
133:       } else {
134:         PetscViewerASCIIPrintf(viewer,"  Inverse square root: %s",str);
135:         PetscViewerASCIIUseTabs(viewer,PETSC_FALSE);
136:         SlepcSNPrintfScalar(str,50,fn->alpha,PETSC_TRUE);
137:         PetscViewerASCIIPrintf(viewer,"*(%s*x)^(-1/2)\n",str);
138:         PetscViewerASCIIUseTabs(viewer,PETSC_TRUE);
139:       }
140:     }
141:   }
142:   return(0);
143: }

147: PETSC_EXTERN PetscErrorCode FNCreate_Invsqrt(FN fn)
148: {
150:   fn->ops->evaluatefunction       = FNEvaluateFunction_Invsqrt;
151:   fn->ops->evaluatederivative     = FNEvaluateDerivative_Invsqrt;
152:   fn->ops->evaluatefunctionmat    = FNEvaluateFunctionMat_Invsqrt;
153:   fn->ops->evaluatefunctionmatvec = FNEvaluateFunctionMatVec_Invsqrt;
154:   fn->ops->view                   = FNView_Invsqrt;
155:   return(0);
156: }