Actual source code: fncombine.c

slepc-3.7.2 2016-07-19
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-2016, 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:       if (b==0.0) SETERRQ(PETSC_COMM_SELF,1,"Function not defined in the requested value");
 60:       *y = a/b;
 61:       break;
 62:     case FN_COMBINE_COMPOSE:
 63:       FNEvaluateFunction(ctx->f2,a,y);
 64:       break;
 65:   }
 66:   return(0);
 67: }

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

 78:   switch (ctx->comb) {
 79:     case FN_COMBINE_ADD:
 80:       FNEvaluateDerivative(ctx->f1,x,&ap);
 81:       FNEvaluateDerivative(ctx->f2,x,&bp);
 82:       *yp = ap+bp;
 83:       break;
 84:     case FN_COMBINE_MULTIPLY:
 85:       FNEvaluateDerivative(ctx->f1,x,&ap);
 86:       FNEvaluateDerivative(ctx->f2,x,&bp);
 87:       FNEvaluateFunction(ctx->f1,x,&a);
 88:       FNEvaluateFunction(ctx->f2,x,&b);
 89:       *yp = ap*b+a*bp;
 90:       break;
 91:     case FN_COMBINE_DIVIDE:
 92:       FNEvaluateDerivative(ctx->f1,x,&ap);
 93:       FNEvaluateDerivative(ctx->f2,x,&bp);
 94:       FNEvaluateFunction(ctx->f1,x,&a);
 95:       FNEvaluateFunction(ctx->f2,x,&b);
 96:       if (b==0.0) SETERRQ(PETSC_COMM_SELF,1,"Derivative not defined in the requested value");
 97:       *yp = (ap*b-a*bp)/(b*b);
 98:       break;
 99:     case FN_COMBINE_COMPOSE:
100:       FNEvaluateFunction(ctx->f1,x,&a);
101:       FNEvaluateDerivative(ctx->f1,x,&ap);
102:       FNEvaluateDerivative(ctx->f2,a,yp);
103:       *yp *= ap;
104:       break;
105:   }
106:   return(0);
107: }

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

125:   FN_AllocateWorkMat(fn,A,&W);
126:   MatDenseGetArray(A,&Aa);
127:   MatDenseGetArray(B,&Ba);
128:   MatDenseGetArray(W,&Wa);
129:   MatGetSize(A,&m,NULL);
130:   PetscBLASIntCast(m,&n);
131:   ld  = n;
132:   ld2 = ld*ld;

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

163:   MatDenseRestoreArray(A,&Aa);
164:   MatDenseRestoreArray(B,&Ba);
165:   MatDenseRestoreArray(W,&Wa);
166:   FN_FreeWorkMat(fn,&W);
167:   return(0);
168: #endif
169: }

173: PetscErrorCode FNEvaluateFunctionMatVec_Combine(FN fn,Mat A,Vec v)
174: {
175: #if defined(PETSC_MISSING_LAPACK_GESV)
177:   SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"GESV - Lapack routines are unavailable");
178: #else
180:   FN_COMBINE     *ctx = (FN_COMBINE*)fn->data;
181:   PetscScalar    *va,*Za;
182:   PetscBLASInt   n,ld,*ipiv,info,one=1;
183:   PetscInt       m;
184:   Mat            Z;
185:   Vec            w;

188:   MatGetSize(A,&m,NULL);
189:   PetscBLASIntCast(m,&n);
190:   ld = n;

192:   switch (ctx->comb) {
193:     case FN_COMBINE_ADD:
194:       VecDuplicate(v,&w);
195:       FNEvaluateFunctionMatVec(ctx->f1,A,w);
196:       FNEvaluateFunctionMatVec(ctx->f2,A,v);
197:       VecAXPY(v,1.0,w);
198:       VecDestroy(&w);
199:       break;
200:     case FN_COMBINE_MULTIPLY:
201:       VecDuplicate(v,&w);
202:       FN_AllocateWorkMat(fn,A,&Z);
203:       FNEvaluateFunctionMat(ctx->f1,A,Z);
204:       FNEvaluateFunctionMatVec(ctx->f2,A,w);
205:       MatMult(Z,w,v);
206:       FN_FreeWorkMat(fn,&Z);
207:       VecDestroy(&w);
208:       break;
209:     case FN_COMBINE_DIVIDE:
210:       VecDuplicate(v,&w);
211:       FN_AllocateWorkMat(fn,A,&Z);
212:       FNEvaluateFunctionMat(ctx->f2,A,Z);
213:       FNEvaluateFunctionMatVec(ctx->f1,A,v);
214:       PetscMalloc1(ld,&ipiv);
215:       MatDenseGetArray(Z,&Za);
216:       VecGetArray(v,&va);
217:       PetscStackCallBLAS("LAPACKgesv",LAPACKgesv_(&n,&one,Za,&ld,ipiv,va,&ld,&info));
218:       if (info) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in Lapack xGESV %d",info);
219:       VecRestoreArray(v,&va);
220:       MatDenseRestoreArray(Z,&Za);
221:       PetscFree(ipiv);
222:       FN_FreeWorkMat(fn,&Z);
223:       VecDestroy(&w);
224:       break;
225:     case FN_COMBINE_COMPOSE:
226:       FN_AllocateWorkMat(fn,A,&Z);
227:       FNEvaluateFunctionMat(ctx->f1,A,Z);
228:       FNEvaluateFunctionMatVec(ctx->f2,Z,v);
229:       FN_FreeWorkMat(fn,&Z);
230:       break;
231:   }
232:   return(0);
233: #endif
234: }

238: PetscErrorCode FNView_Combine(FN fn,PetscViewer viewer)
239: {
241:   FN_COMBINE     *ctx = (FN_COMBINE*)fn->data;
242:   PetscBool      isascii;

245:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);
246:   if (isascii) {
247:     switch (ctx->comb) {
248:       case FN_COMBINE_ADD:
249:         PetscViewerASCIIPrintf(viewer,"  Two added functions f1+f2\n");
250:         break;
251:       case FN_COMBINE_MULTIPLY:
252:         PetscViewerASCIIPrintf(viewer,"  Two multiplied functions f1*f2\n");
253:         break;
254:       case FN_COMBINE_DIVIDE:
255:         PetscViewerASCIIPrintf(viewer,"  A quotient of two functions f1/f2\n");
256:         break;
257:       case FN_COMBINE_COMPOSE:
258:         PetscViewerASCIIPrintf(viewer,"  Two composed functions f2(f1(.))\n");
259:         break;
260:     }
261:     PetscViewerASCIIPushTab(viewer);
262:     FNView(ctx->f1,viewer);
263:     FNView(ctx->f2,viewer);
264:     PetscViewerASCIIPopTab(viewer);
265:   }
266:   return(0);
267: }

271: static PetscErrorCode FNCombineSetChildren_Combine(FN fn,FNCombineType comb,FN f1,FN f2)
272: {
274:   FN_COMBINE     *ctx = (FN_COMBINE*)fn->data;

277:   ctx->comb = comb;
278:   PetscObjectReference((PetscObject)f1);
279:   FNDestroy(&ctx->f1);
280:   ctx->f1 = f1;
281:   PetscLogObjectParent((PetscObject)fn,(PetscObject)ctx->f1);
282:   PetscObjectReference((PetscObject)f2);
283:   FNDestroy(&ctx->f2);
284:   ctx->f2 = f2;
285:   PetscLogObjectParent((PetscObject)fn,(PetscObject)ctx->f2);
286:   return(0);
287: }

291: /*@
292:    FNCombineSetChildren - Sets the two child functions that constitute this
293:    combined function, and the way they must be combined.

295:    Logically Collective on FN

297:    Input Parameters:
298: +  fn   - the math function context
299: .  comb - how to combine the functions (addition, multiplication, division or composition)
300: .  f1   - first function
301: -  f2   - second function

303:    Level: intermediate

305: .seealso: FNCombineGetChildren()
306: @*/
307: PetscErrorCode FNCombineSetChildren(FN fn,FNCombineType comb,FN f1,FN f2)
308: {

316:   PetscTryMethod(fn,"FNCombineSetChildren_C",(FN,FNCombineType,FN,FN),(fn,comb,f1,f2));
317:   return(0);
318: }

322: static PetscErrorCode FNCombineGetChildren_Combine(FN fn,FNCombineType *comb,FN *f1,FN *f2)
323: {
325:   FN_COMBINE     *ctx = (FN_COMBINE*)fn->data;

328:   if (comb) *comb = ctx->comb;
329:   if (f1) {
330:     if (!ctx->f1) {
331:       FNCreate(PetscObjectComm((PetscObject)fn),&ctx->f1);
332:       PetscLogObjectParent((PetscObject)fn,(PetscObject)ctx->f1);
333:     }
334:     *f1 = ctx->f1;
335:   }
336:   if (f2) {
337:     if (!ctx->f2) {
338:       FNCreate(PetscObjectComm((PetscObject)fn),&ctx->f2);
339:       PetscLogObjectParent((PetscObject)fn,(PetscObject)ctx->f2);
340:     }
341:     *f2 = ctx->f2;
342:   }
343:   return(0);
344: }

348: /*@
349:    FNCombineGetChildren - Gets the two child functions that constitute this
350:    combined function, and the way they are combined.

352:    Not Collective

354:    Input Parameter:
355: .  fn   - the math function context

357:    Output Parameters:
358: +  comb - how to combine the functions (addition, multiplication, division or composition)
359: .  f1   - first function
360: -  f2   - second function

362:    Level: intermediate

364: .seealso: FNCombineSetChildren()
365: @*/
366: PetscErrorCode FNCombineGetChildren(FN fn,FNCombineType *comb,FN *f1,FN *f2)
367: {

372:   PetscUseMethod(fn,"FNCombineGetChildren_C",(FN,FNCombineType*,FN*,FN*),(fn,comb,f1,f2));
373:   return(0);
374: }

378: PetscErrorCode FNDuplicate_Combine(FN fn,MPI_Comm comm,FN *newfn)
379: {
381:   FN_COMBINE     *ctx = (FN_COMBINE*)fn->data,*ctx2;

384:   PetscNewLog(*newfn,&ctx2);
385:   (*newfn)->data = (void*)ctx2;
386:   ctx2->comb = ctx->comb;
387:   FNDuplicate(ctx->f1,comm,&ctx2->f1);
388:   FNDuplicate(ctx->f2,comm,&ctx2->f2);
389:   return(0);
390: }

394: PetscErrorCode FNDestroy_Combine(FN fn)
395: {
397:   FN_COMBINE     *ctx = (FN_COMBINE*)fn->data;

400:   FNDestroy(&ctx->f1);
401:   FNDestroy(&ctx->f2);
402:   PetscFree(fn->data);
403:   PetscObjectComposeFunction((PetscObject)fn,"FNCombineSetChildren_C",NULL);
404:   PetscObjectComposeFunction((PetscObject)fn,"FNCombineGetChildren_C",NULL);
405:   return(0);
406: }

410: PETSC_EXTERN PetscErrorCode FNCreate_Combine(FN fn)
411: {
413:   FN_COMBINE     *ctx;

416:   PetscNewLog(fn,&ctx);
417:   fn->data = (void*)ctx;

419:   fn->ops->evaluatefunction       = FNEvaluateFunction_Combine;
420:   fn->ops->evaluatederivative     = FNEvaluateDerivative_Combine;
421:   fn->ops->evaluatefunctionmat    = FNEvaluateFunctionMat_Combine;
422:   fn->ops->evaluatefunctionmatvec = FNEvaluateFunctionMatVec_Combine;
423:   fn->ops->view                   = FNView_Combine;
424:   fn->ops->duplicate              = FNDuplicate_Combine;
425:   fn->ops->destroy                = FNDestroy_Combine;
426:   PetscObjectComposeFunction((PetscObject)fn,"FNCombineSetChildren_C",FNCombineSetChildren_Combine);
427:   PetscObjectComposeFunction((PetscObject)fn,"FNCombineGetChildren_C",FNCombineGetChildren_Combine);
428:   return(0);
429: }