Actual source code: fncombine.c

slepc-3.6.1 2015-09-03
Report Typos and Errors
  1: /*
  2:    A function that is obtained by combining two other functions (either by
  3:    addition, multiplication, division or composition)

  5:       addition:          f(x) = f1(x)+f2(x)
  6:       multiplication:    f(x) = f1(x)*f2(x)
  7:       division:          f(x) = f1(x)/f2(x)      f(A) = f2(A)\f1(A)
  8:       composition:       f(x) = f2(f1(x))

 10:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 11:    SLEPc - Scalable Library for Eigenvalue Problem Computations
 12:    Copyright (c) 2002-2015, Universitat Politecnica de Valencia, Spain

 14:    This file is part of SLEPc.

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

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

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

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

 33: typedef struct {
 34:   FN            f1,f2;    /* functions */
 35:   FNCombineType comb;     /* how the functions are combined */
 36: } FN_COMBINE;

 40: PetscErrorCode FNEvaluateFunction_Combine(FN fn,PetscScalar x,PetscScalar *y)
 41: {
 43:   FN_COMBINE     *ctx = (FN_COMBINE*)fn->data;
 44:   PetscScalar    a,b;

 47:   FNEvaluateFunction(ctx->f1,x,&a);
 48:   switch (ctx->comb) {
 49:     case FN_COMBINE_ADD:
 50:       FNEvaluateFunction(ctx->f2,x,&b);
 51:       *y = a+b;
 52:       break;
 53:     case FN_COMBINE_MULTIPLY:
 54:       FNEvaluateFunction(ctx->f2,x,&b);
 55:       *y = a*b;
 56:       break;
 57:     case FN_COMBINE_DIVIDE:
 58:       FNEvaluateFunction(ctx->f2,x,&b);
 59:       *y = a/b;
 60:       break;
 61:     case FN_COMBINE_COMPOSE:
 62:       FNEvaluateFunction(ctx->f2,a,y);
 63:       break;
 64:   }
 65:   return(0);
 66: }

 70: PetscErrorCode FNEvaluateDerivative_Combine(FN fn,PetscScalar x,PetscScalar *yp)
 71: {
 73:   FN_COMBINE     *ctx = (FN_COMBINE*)fn->data;
 74:   PetscScalar    a,b,ap,bp;

 77:   switch (ctx->comb) {
 78:     case FN_COMBINE_ADD:
 79:       FNEvaluateDerivative(ctx->f1,x,&ap);
 80:       FNEvaluateDerivative(ctx->f2,x,&bp);
 81:       *yp = ap+bp;
 82:       break;
 83:     case FN_COMBINE_MULTIPLY:
 84:       FNEvaluateDerivative(ctx->f1,x,&ap);
 85:       FNEvaluateDerivative(ctx->f2,x,&bp);
 86:       FNEvaluateFunction(ctx->f1,x,&a);
 87:       FNEvaluateFunction(ctx->f2,x,&b);
 88:       *yp = ap*b+a*bp;
 89:       break;
 90:     case FN_COMBINE_DIVIDE:
 91:       FNEvaluateDerivative(ctx->f1,x,&ap);
 92:       FNEvaluateDerivative(ctx->f2,x,&bp);
 93:       FNEvaluateFunction(ctx->f1,x,&a);
 94:       FNEvaluateFunction(ctx->f2,x,&b);
 95:       *yp = (ap*b-a*bp)/(b*b);
 96:       break;
 97:     case FN_COMBINE_COMPOSE:
 98:       FNEvaluateFunction(ctx->f1,x,&a);
 99:       FNEvaluateDerivative(ctx->f1,x,&ap);
100:       FNEvaluateDerivative(ctx->f2,a,yp);
101:       *yp *= ap;
102:       break;
103:   }
104:   return(0);
105: }

109: PetscErrorCode FNEvaluateFunctionMat_Combine(FN fn,Mat A,Mat B)
110: {
111: #if defined(PETSC_MISSING_LAPACK_GESV)
113:   SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"GESV - Lapack routines are unavailable");
114: #else
116:   FN_COMBINE     *ctx = (FN_COMBINE*)fn->data;
117:   PetscScalar    *Aa,*Ba,*Wa,*Za,one=1.0,zero=0.0;
118:   PetscBLASInt   n,ld,ld2,inc=1,*ipiv,info;
119:   PetscInt       m;
120:   Mat            W,Z;

123:   MatDuplicate(A,MAT_DO_NOT_COPY_VALUES,&W);
124:   MatDenseGetArray(A,&Aa);
125:   MatDenseGetArray(B,&Ba);
126:   MatDenseGetArray(W,&Wa);
127:   MatGetSize(A,&m,NULL);
128:   PetscBLASIntCast(m,&n);
129:   ld  = n;
130:   ld2 = ld*ld;

132:   switch (ctx->comb) {
133:     case FN_COMBINE_ADD:
134:       FNEvaluateFunctionMat(ctx->f1,A,W);
135:       FNEvaluateFunctionMat(ctx->f2,A,B);
136:       PetscStackCallBLAS("BLASaxpy",BLASaxpy_(&ld2,&one,Wa,&inc,Ba,&inc));
137:       break;
138:     case FN_COMBINE_MULTIPLY:
139:       MatDuplicate(A,MAT_DO_NOT_COPY_VALUES,&Z);
140:       MatDenseGetArray(Z,&Za);
141:       FNEvaluateFunctionMat(ctx->f1,A,W);
142:       FNEvaluateFunctionMat(ctx->f2,A,Z);
143:       PetscStackCallBLAS("BLASgemm",BLASgemm_("N","N",&n,&n,&n,&one,Wa,&ld,Za,&ld,&zero,Ba,&ld));
144:       MatDenseRestoreArray(Z,&Za);
145:       MatDestroy(&Z);
146:       break;
147:     case FN_COMBINE_DIVIDE:
148:       FNEvaluateFunctionMat(ctx->f1,A,B);
149:       FNEvaluateFunctionMat(ctx->f2,A,W);
150:       PetscMalloc1(ld,&ipiv);
151:       PetscStackCallBLAS("LAPACKgesv",LAPACKgesv_(&n,&n,Wa,&ld,ipiv,Ba,&ld,&info));
152:       if (info) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in Lapack xGESV %d",info);
153:       PetscFree(ipiv);
154:       break;
155:     case FN_COMBINE_COMPOSE:
156:       FNEvaluateFunctionMat(ctx->f1,A,W);
157:       FNEvaluateFunctionMat(ctx->f2,W,B);
158:       break;
159:   }

161:   MatDenseRestoreArray(A,&Aa);
162:   MatDenseRestoreArray(B,&Ba);
163:   MatDenseRestoreArray(W,&Wa);
164:   MatDestroy(&W);
165:   return(0);
166: #endif
167: }

171: PetscErrorCode FNView_Combine(FN fn,PetscViewer viewer)
172: {
174:   FN_COMBINE     *ctx = (FN_COMBINE*)fn->data;
175:   PetscBool      isascii;

178:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);
179:   if (isascii) {
180:     switch (ctx->comb) {
181:       case FN_COMBINE_ADD:
182:         PetscViewerASCIIPrintf(viewer,"  Two added functions f1+f2\n");
183:         break;
184:       case FN_COMBINE_MULTIPLY:
185:         PetscViewerASCIIPrintf(viewer,"  Two multiplied functions f1*f2\n");
186:         break;
187:       case FN_COMBINE_DIVIDE:
188:         PetscViewerASCIIPrintf(viewer,"  A quotient of two functions f1/f2\n");
189:         break;
190:       case FN_COMBINE_COMPOSE:
191:         PetscViewerASCIIPrintf(viewer,"  Two composed functions f2(f1(.))\n");
192:         break;
193:     }
194:     PetscViewerASCIIPushTab(viewer);
195:     FNView(ctx->f1,viewer);
196:     FNView(ctx->f2,viewer);
197:     PetscViewerASCIIPopTab(viewer);
198:   }
199:   return(0);
200: }

204: static PetscErrorCode FNCombineSetChildren_Combine(FN fn,FNCombineType comb,FN f1,FN f2)
205: {
207:   FN_COMBINE     *ctx = (FN_COMBINE*)fn->data;

210:   ctx->comb = comb;
211:   PetscObjectReference((PetscObject)f1);
212:   FNDestroy(&ctx->f1);
213:   ctx->f1 = f1;
214:   PetscLogObjectParent((PetscObject)fn,(PetscObject)ctx->f1);
215:   PetscObjectReference((PetscObject)f2);
216:   FNDestroy(&ctx->f2);
217:   ctx->f2 = f2;
218:   PetscLogObjectParent((PetscObject)fn,(PetscObject)ctx->f2);
219:   return(0);
220: }

224: /*@
225:    FNCombineSetChildren - Sets the two child functions that constitute this
226:    combined function, and the way they must be combined.

228:    Logically Collective on FN

230:    Input Parameters:
231: +  fn   - the math function context
232: .  comb - how to combine the functions (addition, multiplication, division or composition)
233: .  f1   - first function
234: -  f2   - second function

236:    Level: intermediate

238: .seealso: FNCombineGetChildren()
239: @*/
240: PetscErrorCode FNCombineSetChildren(FN fn,FNCombineType comb,FN f1,FN f2)
241: {

249:   PetscTryMethod(fn,"FNCombineSetChildren_C",(FN,FNCombineType,FN,FN),(fn,comb,f1,f2));
250:   return(0);
251: }

255: static PetscErrorCode FNCombineGetChildren_Combine(FN fn,FNCombineType *comb,FN *f1,FN *f2)
256: {
258:   FN_COMBINE     *ctx = (FN_COMBINE*)fn->data;

261:   if (comb) *comb = ctx->comb;
262:   if (f1) {
263:     if (!ctx->f1) {
264:       FNCreate(PetscObjectComm((PetscObject)fn),&ctx->f1);
265:       PetscLogObjectParent((PetscObject)fn,(PetscObject)ctx->f1);
266:     }
267:     *f1 = ctx->f1;
268:   }
269:   if (f2) {
270:     if (!ctx->f2) {
271:       FNCreate(PetscObjectComm((PetscObject)fn),&ctx->f2);
272:       PetscLogObjectParent((PetscObject)fn,(PetscObject)ctx->f2);
273:     }
274:     *f2 = ctx->f2;
275:   }
276:   return(0);
277: }

281: /*@
282:    FNCombineGetChildren - Gets the two child functions that constitute this
283:    combined function, and the way they are combined.

285:    Not Collective

287:    Input Parameter:
288: .  fn   - the math function context

290:    Output Parameters:
291: -  comb - how to combine the functions (addition, multiplication, division or composition)
292: .  f1   - first function
293: -  f2   - second function

295:    Level: intermediate

297: .seealso: FNCombineSetChildren()
298: @*/
299: PetscErrorCode FNCombineGetChildren(FN fn,FNCombineType *comb,FN *f1,FN *f2)
300: {

305:   PetscTryMethod(fn,"FNCombineGetChildren_C",(FN,FNCombineType*,FN*,FN*),(fn,comb,f1,f2));
306:   return(0);
307: }

311: PetscErrorCode FNDuplicate_Combine(FN fn,MPI_Comm comm,FN *newfn)
312: {
314:   FN_COMBINE     *ctx = (FN_COMBINE*)fn->data,*ctx2;

317:   PetscNewLog(*newfn,&ctx2);
318:   (*newfn)->data = (void*)ctx2;
319:   ctx2->comb = ctx->comb;
320:   FNDuplicate(ctx->f1,comm,&ctx2->f1);
321:   FNDuplicate(ctx->f2,comm,&ctx2->f2);
322:   return(0);
323: }

327: PetscErrorCode FNDestroy_Combine(FN fn)
328: {
330:   FN_COMBINE     *ctx = (FN_COMBINE*)fn->data;

333:   FNDestroy(&ctx->f1);
334:   FNDestroy(&ctx->f2);
335:   PetscFree(fn->data);
336:   PetscObjectComposeFunction((PetscObject)fn,"FNCombineSetChildren_C",NULL);
337:   PetscObjectComposeFunction((PetscObject)fn,"FNCombineGetChildren_C",NULL);
338:   return(0);
339: }

343: PETSC_EXTERN PetscErrorCode FNCreate_Combine(FN fn)
344: {
346:   FN_COMBINE     *ctx;

349:   PetscNewLog(fn,&ctx);
350:   fn->data = (void*)ctx;

352:   fn->ops->evaluatefunction    = FNEvaluateFunction_Combine;
353:   fn->ops->evaluatederivative  = FNEvaluateDerivative_Combine;
354:   fn->ops->evaluatefunctionmat = FNEvaluateFunctionMat_Combine;
355:   fn->ops->view                = FNView_Combine;
356:   fn->ops->duplicate           = FNDuplicate_Combine;
357:   fn->ops->destroy             = FNDestroy_Combine;
358:   PetscObjectComposeFunction((PetscObject)fn,"FNCombineSetChildren_C",FNCombineSetChildren_Combine);
359:   PetscObjectComposeFunction((PetscObject)fn,"FNCombineGetChildren_C",FNCombineGetChildren_Combine);
360:   return(0);
361: }