1 #ifdef PETSC_RCS_HEADER 2 static char vcid[] = "$Id: snesmfj.c,v 1.71 1998/10/31 23:58:40 bsmith Exp bsmith $"; 3 #endif 4 5 #include "src/snes/snesimpl.h" /*I "snes.h" I*/ 6 #include "pinclude/pviewer.h" 7 8 typedef struct _p_MatSNESFDMF_HCtx *MatSNESFDMF_HCtx; 9 typedef struct _p_MatSNESFDMFCtx *MatSNESFDMFCtx; 10 11 typedef struct { 12 int (*compute)(MatSNESFDMFCtx,Vec,Vec,Scalar *); 13 int (*view)(MatSNESFDMFCtx); 14 int (*destroy)(MatSNESFDMFCtx); 15 } MFHOps; 16 17 struct _p_MatSNESFDMF_HCtx { 18 MFHOps *ops; 19 void *data; 20 }; 21 22 struct _p_MatSNESFDMFCtx { /* context for default matrix-free SNES */ 23 SNES snes; /* SNES context */ 24 Vec w; /* work vector */ 25 PCNullSpace sp; /* null space context */ 26 double error_rel; /* square root of relative error in computing function */ 27 double umin; /* minimum allowable u'a value relative to |u|_1 */ 28 Scalar currenth; /* last differencing parameter used */ 29 Scalar *historyh; /* history of h */ 30 int ncurrenth,maxcurrenth; 31 MatSNESFDMF_HCtx hctx; 32 }; 33 34 #undef __FUNC__ 35 #define __FUNC__ "MatSNESFDMFSetType" 36 int MatSNESFDMFDestroy_Private(Mat mat) 37 { 38 int ierr; 39 MatSNESFDMFCtx ctx; 40 41 PetscFunctionBegin; 42 ierr = MatShellGetContext(mat,(void **)&ctx);CHKERRQ(ierr); 43 44 45 #undef __FUNC__ 46 #define __FUNC__ "MatSNESFDMFDestroy_Private" 47 int MatSNESFDMFDestroy_Private(Mat mat) 48 { 49 int ierr; 50 MatSNESFDMFCtx ctx; 51 52 PetscFunctionBegin; 53 ierr = MatShellGetContext(mat,(void **)&ctx);CHKERRQ(ierr); 54 ierr = VecDestroy(ctx->w); CHKERRQ(ierr); 55 if (ctx->hctx->ops->destroy) {ierr = (*ctx->hctx->ops->destroy)(ctx);CHKERRQ(ierr);} 56 if (ctx->sp) {ierr = PCNullSpaceDestroy(ctx->sp);CHKERRQ(ierr);} 57 PetscFree(ctx); 58 PetscFunctionReturn(0); 59 } 60 61 #undef __FUNC__ 62 #define __FUNC__ "MatSNESFDMFView_Private" 63 /* 64 MatSNESFDMFView_Private - Views matrix-free parameters. 65 66 */ 67 int MatSNESFDMFView_Private(Mat J,Viewer viewer) 68 { 69 int ierr; 70 MatSNESFDMFCtx ctx; 71 MPI_Comm comm; 72 FILE *fd; 73 ViewerType vtype; 74 75 PetscFunctionBegin; 76 PetscObjectGetComm((PetscObject)J,&comm); 77 ierr = MatShellGetContext(J,(void **)&ctx); CHKERRQ(ierr); 78 ierr = ViewerGetType(viewer,&vtype); CHKERRQ(ierr); 79 ierr = ViewerASCIIGetPointer(viewer,&fd); CHKERRQ(ierr); 80 if (vtype == ASCII_FILE_VIEWER || vtype == ASCII_FILES_VIEWER) { 81 PetscFPrintf(comm,fd," SNES matrix-free approximation:\n"); 82 PetscFPrintf(comm,fd," err=%g (relative error in function evaluation)\n",ctx->error_rel); 83 PetscFPrintf(comm,fd," umin=%g (minimum iterate parameter)\n",ctx->umin); 84 } else { 85 SETERRQ(1,1,"Viewer type not supported for this object"); 86 } 87 PetscFunctionReturn(0); 88 } 89 90 91 #undef __FUNC__ 92 #define __FUNC__ "MatSNESFDMFMult_Private" 93 /* 94 MatSNESFDMFMult_Private - Default matrix-free form for Jacobian-vector 95 product, y = F'(u)*a: 96 97 y = ( F(u + ha) - F(u) ) /h, 98 where F = nonlinear function, as set by SNESSetFunction() 99 u = current iterate 100 h = difference interval 101 */ 102 int MatSNESFDMFMult_Private(Mat mat,Vec a,Vec y) 103 { 104 MatSNESFDMFCtx ctx; 105 SNES snes; 106 Scalar h, mone = -1.0; 107 Vec w,U,F; 108 int ierr, (*eval_fct)(SNES,Vec,Vec); 109 110 PetscFunctionBegin; 111 PLogEventBegin(MAT_MatrixFreeMult,a,y,0,0); 112 113 ierr = MatShellGetContext(mat,(void **)&ctx); CHKERRQ(ierr); 114 snes = ctx->snes; 115 w = ctx->w; 116 117 /* We log matrix-free matrix-vector products separately, so that we can 118 separate the performance monitoring from the cases that use conventional 119 storage. We may eventually modify event logging to associate events 120 with particular objects, hence alleviating the more general problem. */ 121 122 ierr = SNESGetSolution(snes,&U); CHKERRQ(ierr); 123 if (snes->method_class == SNES_NONLINEAR_EQUATIONS) { 124 eval_fct = SNESComputeFunction; 125 ierr = SNESGetFunction(snes,&F); CHKERRQ(ierr); 126 } else if (snes->method_class == SNES_UNCONSTRAINED_MINIMIZATION) { 127 eval_fct = SNESComputeGradient; 128 ierr = SNESGetGradient(snes,&F); CHKERRQ(ierr); 129 } else SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Invalid method class"); 130 131 ierr = (*ctx->hctx->ops->compute)(ctx,U,a,&h); CHKERRQ(ierr); 132 133 /* keep a record of the current differencing parameter h */ 134 ctx->currenth = h; 135 #if defined(USE_PETSC_COMPLEX) 136 PLogInfo(mat,"Current differencing parameter: %g + %g i\n",PetscReal(h),PetscImaginary(h)); 137 #else 138 PLogInfo(mat,"Current differencing parameter: %g\n",h); 139 #endif 140 if (ctx->historyh && ctx->ncurrenth < ctx->maxcurrenth) { 141 ctx->historyh[ctx->ncurrenth++] = h; 142 } 143 144 /* Evaluate function at F(u + ha) */ 145 ierr = VecWAXPY(&h,a,U,w); CHKERRQ(ierr); 146 ierr = eval_fct(snes,w,y); CHKERRQ(ierr); 147 148 ierr = VecAXPY(&mone,F,y); CHKERRQ(ierr); 149 h = 1.0/h; 150 ierr = VecScale(&h,y); CHKERRQ(ierr); 151 if (ctx->sp) {ierr = PCNullSpaceRemove(ctx->sp,y); CHKERRQ(ierr);} 152 153 PLogEventEnd(MAT_MatrixFreeMult,a,y,0,0); 154 PetscFunctionReturn(0); 155 } 156 157 extern int MatFDMFDefaultComputeH(SNESFDMFCtx,Vec,Vec,Scalar *); 158 159 #undef __FUNC__ 160 #define __FUNC__ "MatCreateSNESFDMF" 161 /*@C 162 MatCreateSNESFDMF - Creates a matrix-free matrix 163 context for use with a SNES solver. This matrix can be used as 164 the Jacobian argument for the routine SNESSetJacobian(). 165 166 Collective on SNES and Vec 167 168 Input Parameters: 169 + snes - the SNES context 170 - x - vector where SNES solution is to be stored. 171 172 Output Parameter: 173 . J - the matrix-free matrix 174 175 Notes: 176 The matrix-free matrix context merely contains the function pointers 177 and work space for performing finite difference approximations of 178 Jacobian-vector products, J(u)*a, via 179 180 .vb 181 J(u)*a = [J(u+h*a) - J(u)]/h where 182 h = error_rel*u'a/||a||^2 if |u'a| > umin*||a||_{1} 183 = error_rel*umin*sign(u'a)*||a||_{1}/||a||^2 otherwise 184 where 185 error_rel = square root of relative error in function evaluation 186 umin = minimum iterate parameter 187 .ve 188 189 The user can set these parameters via MatSNESFDMFSetParameters(). 190 See the nonlinear solvers chapter of the users manual for details. 191 192 The user should call MatDestroy() when finished with the matrix-free 193 matrix context. 194 195 Options Database Keys: 196 + -snes_mf_err <error_rel> - Sets error_rel 197 . -snes_mf_unim <umin> - Sets umin 198 - -snes_mf_ksp_monitor - KSP monitor routine that prints differencing h 199 200 .keywords: SNES, default, matrix-free, create, matrix 201 202 .seealso: MatDestroy(), MatSNESFDMFSetParameters(), 203 MatSNESFDMFSetHHistory(), MatSNESFDMFResetHHistory(), 204 MatSNESFDMFGetH(),MatSNESFDMFKSPMonitor() 205 206 @*/ 207 int MatCreateSNESFDMF(SNES snes,Vec x, Mat *J) 208 { 209 MPI_Comm comm; 210 MatSNESFDMFCtx mfctx; 211 int n, nloc, ierr, flg; 212 char p[64]; 213 214 PetscFunctionBegin; 215 mfctx = (MatSNESFDMFCtx) PetscMalloc(sizeof(struct _p_MatSNESFDMFCtx)); CHKPTRQ(mfctx); 216 PLogObjectMemory(snes,sizeof(MatSNESFDMFCtx)); 217 mfctx->sp = 0; 218 mfctx->snes = snes; 219 mfctx->error_rel = 1.e-8; /* assumes double precision */ 220 mfctx->umin = 1.e-6; 221 mfctx->currenth = 0.0; 222 mfctx->historyh = PETSC_NULL; 223 mfctx->ncurrenth = 0; 224 mfctx->maxcurrenth = 0; 225 mfctx->hctx = 0; 226 mfctx->hctx->ops->compute = MatFDMFDefaultComputeH; 227 mfctx->hctx->ops->destroy = 0; 228 mfctx->hctx->ops->view = 0; 229 230 ierr = OptionsGetDouble(snes->prefix,"-snes_mf_err",&mfctx->error_rel,&flg); CHKERRQ(ierr); 231 ierr = OptionsGetDouble(snes->prefix,"-snes_mf_umin",&mfctx->umin,&flg); CHKERRQ(ierr); 232 ierr = OptionsHasName(PETSC_NULL,"-help",&flg); CHKERRQ(ierr); 233 PetscStrcpy(p,"-"); 234 if (snes->prefix) PetscStrcat(p,snes->prefix); 235 if (flg) { 236 (*PetscHelpPrintf)(snes->comm," %ssnes_mf_err <err>: set sqrt rel error in function (default %g)\n",p,mfctx->error_rel); 237 (*PetscHelpPrintf)(snes->comm," %ssnes_mf_umin <umin> see users manual (default %g)\n",p,mfctx->umin); 238 } 239 ierr = VecDuplicate(x,&mfctx->w); CHKERRQ(ierr); 240 ierr = PetscObjectGetComm((PetscObject)x,&comm); CHKERRQ(ierr); 241 ierr = VecGetSize(x,&n); CHKERRQ(ierr); 242 ierr = VecGetLocalSize(x,&nloc); CHKERRQ(ierr); 243 ierr = MatCreateShell(comm,nloc,n,n,n,mfctx,J); CHKERRQ(ierr); 244 ierr = MatShellSetOperation(*J,MATOP_MULT,(void*)MatSNESFDMFMult_Private);CHKERRQ(ierr); 245 ierr = MatShellSetOperation(*J,MATOP_DESTROY,(void *)MatSNESFDMFDestroy_Private);CHKERRQ(ierr); 246 ierr = MatShellSetOperation(*J,MATOP_VIEW,(void *)MatSNESFDMFView_Private); CHKERRQ(ierr); 247 PLogObjectParent(*J,mfctx->w); 248 PLogObjectParent(snes,*J); 249 PetscFunctionReturn(0); 250 } 251 252 #undef __FUNC__ 253 #define __FUNC__ "MatSNESFDMFGetH" 254 /*@ 255 MatSNESFDMFGetH - Gets the last h that was used as the differencing 256 parameter. 257 258 Not Collective 259 260 Input Parameters: 261 . mat - the matrix obtained with MatCreateSNESFDMF() 262 263 Output Paramter: 264 . h - the differencing step size 265 266 .keywords: SNES, matrix-free, parameters 267 268 .seealso: MatCreateSNESFDMF(),MatSNESFDMFSetHHistory(), 269 MatSNESFDMFResetHHistory(),MatSNESFDMFKSPMonitor() 270 @*/ 271 int MatSNESFDMFGetH(Mat mat,Scalar *h) 272 { 273 MatSNESFDMFCtx ctx; 274 int ierr; 275 276 PetscFunctionBegin; 277 ierr = MatShellGetContext(mat,(void **)&ctx); CHKERRQ(ierr); 278 if (ctx) { 279 *h = ctx->currenth; 280 } 281 PetscFunctionReturn(0); 282 } 283 284 #undef __FUNC__ 285 #define __FUNC__ "MatSNESFDMFKSPMonitor" 286 /* 287 MatSNESFDMFKSPMonitor - A KSP monitor for use with the default PETSc 288 SNES matrix free routines. Prints the h differencing parameter used at each 289 timestep. 290 291 */ 292 int MatSNESFDMFKSPMonitor(KSP ksp,int n,double rnorm,void *dummy) 293 { 294 PC pc; 295 MatSNESFDMFCtx ctx; 296 int ierr; 297 Mat mat; 298 MPI_Comm comm; 299 PetscTruth nonzeroinitialguess; 300 301 PetscFunctionBegin; 302 ierr = PetscObjectGetComm((PetscObject)ksp,&comm);CHKERRQ(ierr); 303 ierr = KSPGetPC(ksp,&pc); CHKERRQ(ierr); 304 ierr = KSPGetInitialGuessNonzero(ksp,&nonzeroinitialguess);CHKERRQ(ierr); 305 ierr = PCGetOperators(pc,&mat,PETSC_NULL,PETSC_NULL); CHKERRQ(ierr); 306 ierr = MatShellGetContext(mat,(void **)&ctx); CHKERRQ(ierr); 307 if (n > 0 || nonzeroinitialguess) { 308 #if defined(USE_PETSC_COMPLEX) 309 PetscPrintf(comm,"%d KSP Residual norm %14.12e h %g + %g i\n",n,rnorm, 310 PetscReal(ctx->currenth),PetscImaginary(ctx->currenth)); 311 #else 312 PetscPrintf(comm,"%d KSP Residual norm %14.12e h %g \n",n,rnorm,ctx->currenth); 313 #endif 314 } else { 315 PetscPrintf(comm,"%d KSP Residual norm %14.12e\n",n,rnorm); 316 } 317 PetscFunctionReturn(0); 318 } 319 320 #undef __FUNC__ 321 #define __FUNC__ "MatSNESFDMFSetParameters" 322 /*@ 323 MatSNESFDMFSetParameters - Sets the parameters for the approximation of 324 matrix-vector products using finite differences. 325 326 Collective on Mat 327 328 Input Parameters: 329 + mat - the matrix free matrix created via MatCreateSNESFDMF() 330 . error_rel - relative error (should be set to the square root of 331 the relative error in the function evaluations) 332 - umin - minimum allowable u-value 333 334 Notes: 335 The default matrix-free matrix-vector product routine computes 336 .vb 337 J(u)*a = [J(u+h*a) - J(u)]/h where 338 h = error_rel*u'a/||a||^2 if |u'a| > umin*||a||_{1} 339 = error_rel*umin*sign(u'a)*||a||_{1}/||a||^2 else 340 .ve 341 342 Options Database Keys: 343 + -snes_mf_err <error_rel> - Sets error_rel 344 - -snes_mf_unim <umin> - Sets umin 345 346 .keywords: SNES, matrix-free, parameters 347 348 .seealso: MatCreateSNESFDMF(),MatSNESFDMFGetH(), 349 MatSNESFDMFSetHHistory(), MatSNESFDMFResetHHistory(), 350 MatSNESFDMFKSPMonitor() 351 @*/ 352 int MatSNESFDMFSetParameters(Mat mat,double error,double umin) 353 { 354 MatSNESFDMFCtx ctx; 355 int ierr; 356 357 PetscFunctionBegin; 358 ierr = MatShellGetContext(mat,(void **)&ctx); CHKERRQ(ierr); 359 if (ctx) { 360 if (error != PETSC_DEFAULT) ctx->error_rel = error; 361 if (umin != PETSC_DEFAULT) ctx->umin = umin; 362 } 363 PetscFunctionReturn(0); 364 } 365 366 #undef __FUNC__ 367 #define __FUNC__ "MatSNESFDMFAddNullSpace" 368 /*@ 369 MatSNESFDMFAddNullSpace - Provides a null space that 370 an operator is supposed to have. Since roundoff will create a 371 small component in the null space, if you know the null space 372 you may have it automatically removed. 373 374 Collective on Mat 375 376 Input Parameters: 377 + J - the matrix-free matrix context 378 . has_cnst - PETSC_TRUE or PETSC_FALSE, indicating if null space has constants 379 . n - number of vectors (excluding constant vector) in null space 380 - vecs - the vectors that span the null space (excluding the constant vector); 381 these vectors must be orthonormal 382 383 .keywords: SNES, matrix-free, null space 384 385 .seealso: MatSNESFDMFGetH(), MatCreateSNESFDMF(), 386 MatSNESFDMFSetHHistory(), MatSNESFDMFResetHHistory(), 387 MatSNESFDMFKSPMonitor(), MatSNESFDMFSetParameters() 388 389 @*/ 390 int MatSNESFDMFAddNullSpace(Mat J,int has_cnst,int n,Vec *vecs) 391 { 392 int ierr; 393 MatSNESFDMFCtx ctx; 394 MPI_Comm comm; 395 396 PetscFunctionBegin; 397 PetscObjectGetComm((PetscObject)J,&comm); 398 399 ierr = MatShellGetContext(J,(void **)&ctx); CHKERRQ(ierr); 400 /* no context indicates that it is not the "matrix free" matrix type */ 401 if (!ctx) PetscFunctionReturn(0); 402 ierr = PCNullSpaceCreate(comm,has_cnst,n,vecs,&ctx->sp); CHKERRQ(ierr); 403 PetscFunctionReturn(0); 404 } 405 406 #undef __FUNC__ 407 #define __FUNC__ "MatSNESFDMFSetHHistory" 408 /*@ 409 MatSNESFDMFSetHHistory - Sets an array to collect a history 410 of the differencing values h computed for the matrix free product 411 412 Collective on Mat 413 414 Input Parameters: 415 + J - the matrix-free matrix context 416 . histroy - space to hold the h history 417 - nhistory - number of entries in history, if more h are generated than 418 nhistory the later ones are discarded 419 420 Notes: 421 Use MatSNESFDMFResetHHistory() to reset the history counter 422 and collect a new batch of h. 423 424 .keywords: SNES, matrix-free, h history, differencing history 425 426 .seealso: MatSNESFDMFGetH(), MatCreateSNESFDMF(), 427 MatSNESFDMFResetHHistory(), 428 MatSNESFDMFKSPMonitor(), MatSNESFDMFSetParameters() 429 430 @*/ 431 int MatSNESFDMFSetHHistory(Mat J,Scalar *history,int nhistory) 432 { 433 int ierr; 434 MatSNESFDMFCtx ctx; 435 436 PetscFunctionBegin; 437 438 ierr = MatShellGetContext(J,(void **)&ctx); CHKERRQ(ierr); 439 /* no context indicates that it is not the "matrix free" matrix type */ 440 if (!ctx) PetscFunctionReturn(0); 441 ctx->historyh = history; 442 ctx->maxcurrenth = nhistory; 443 ctx->currenth = 0; 444 445 PetscFunctionReturn(0); 446 } 447 448 #undef __FUNC__ 449 #define __FUNC__ "MatSNESFDMFResetHHistory" 450 /*@ 451 MatSNESFDMFResetHHistory - Resets the counter to zero to begin 452 collecting a new set of differencing histories. 453 454 Collective on Mat 455 456 Input Parameters: 457 . J - the matrix-free matrix context 458 459 Notes: 460 Use MatSNESFDMFSetHHistory() to create the original history counter 461 462 .keywords: SNES, matrix-free, h history, differencing history 463 464 .seealso: MatSNESFDMFGetH(), MatCreateSNESFDMF(), 465 MatSNESFDMFSetHHistory(), 466 MatSNESFDMFKSPMonitor(), MatSNESFDMFSetParameters() 467 468 @*/ 469 int MatSNESFDMFResetHHistory(Mat J) 470 { 471 int ierr; 472 MatSNESFDMFCtx ctx; 473 474 PetscFunctionBegin; 475 476 ierr = MatShellGetContext(J,(void **)&ctx); CHKERRQ(ierr); 477 /* no context indicates that it is not the "matrix free" matrix type */ 478 if (!ctx) PetscFunctionReturn(0); 479 ctx->currenth = 0; 480 481 PetscFunctionReturn(0); 482 } 483 484 /* ---------------------------------------------------------------------------- */ 485 extern int VecDot_Seq(Vec,Vec,Scalar *); 486 extern int VecNorm_Seq(Vec,NormType,double *); 487 488 typedef struct { 489 double error_rel; /* square root of relative error in computing function */ 490 double umin; /* minimum allowable u'a value relative to |u|_1 */ 491 } MatSNESFDMFDefaultComputeH_struct; 492 493 #undef __FUNC__ 494 #define __FUNC__ "MatSNESFDMFDefaultComputeH" 495 /* 496 MatSNESFDMFDefaultComputeH - Standard PETSc code for 497 computing h with matrix-free finite differences. 498 499 */ 500 int MatFDMFDefaultComputeH(MatSNESFDMFCtx ctx,Vec U,Vec a,Scalar *h) 501 { 502 MPI_Comm comm; 503 double ovalues[3],norm, sum, umin; 504 Scalar dot; 505 int ierr; 506 #if !defined(USE_PETSC_COMPLEX) 507 double values[3]; 508 #endif 509 510 umin = ctx->umin; 511 512 513 /* 514 ierr = VecDot(U,a,&dot); CHKERRQ(ierr); 515 ierr = VecNorm(a,NORM_1,&sum); CHKERRQ(ierr); 516 ierr = VecNorm(a,NORM_2,&norm); CHKERRQ(ierr); 517 */ 518 519 /* 520 Call the Seq Vector routines and then do a single reduction 521 to reduce the number of communications required 522 */ 523 524 #if !defined(USE_PETSC_COMPLEX) 525 PLogEventBegin(VEC_Dot,U,a,0,0); 526 ierr = VecDot_Seq(U,a,ovalues); CHKERRQ(ierr); 527 PLogEventEnd(VEC_Dot,U,a,0,0); 528 PLogEventBegin(VEC_Norm,a,0,0,0); 529 ierr = VecNorm_Seq(a,NORM_1,ovalues+1); CHKERRQ(ierr); 530 ierr = VecNorm_Seq(a,NORM_2,ovalues+2); CHKERRQ(ierr); 531 ovalues[2] = ovalues[2]*ovalues[2]; 532 PLogEventBarrierBegin(VEC_NormBarrier,0,0,0,0,comm); 533 ierr = MPI_Allreduce(ovalues,values,3,MPI_DOUBLE,MPI_SUM,comm );CHKERRQ(ierr); 534 PLogEventBarrierEnd(VEC_NormBarrier,0,0,0,0,comm); 535 dot = values[0]; sum = values[1]; norm = sqrt(values[2]); 536 PLogEventEnd(VEC_Norm,a,0,0,0); 537 #else 538 { 539 Scalar cvalues[3],covalues[3]; 540 541 PLogEventBegin(VEC_Dot,U,a,0,0); 542 ierr = VecDot_Seq(U,a,covalues); CHKERRQ(ierr); 543 PLogEventEnd(VEC_Dot,U,a,0,0); 544 PLogEventBegin(VEC_Norm,a,0,0,0); 545 ierr = VecNorm_Seq(a,NORM_1,ovalues+1); CHKERRQ(ierr); 546 ierr = VecNorm_Seq(a,NORM_2,ovalues+2); CHKERRQ(ierr); 547 covalues[1] = ovalues[1]; 548 covalues[2] = ovalues[2]*ovalues[2]; 549 PLogEventBarrierBegin(VEC_NormBarrier,0,0,0,0,comm); 550 ierr = MPI_Allreduce(covalues,cvalues,6,MPI_DOUBLE,MPI_SUM,comm);CHKERRQ(ierr); 551 PLogEventBarrierBegin(VEC_NormBarrier,0,0,0,0,comm); 552 dot = cvalues[0]; sum = PetscReal(cvalues[1]); norm = sqrt(PetscReal(cvalues[2])); 553 PLogEventEnd(VEC_Norm,a,0,0,0); 554 } 555 #endif 556 557 /* Safeguard for step sizes too small */ 558 if (sum == 0.0) {dot = 1.0; norm = 1.0;} 559 #if defined(USE_PETSC_COMPLEX) 560 else if (PetscAbsScalar(dot) < umin*sum && PetscReal(dot) >= 0.0) dot = umin*sum; 561 else if (PetscAbsScalar(dot) < 0.0 && PetscReal(dot) > -umin*sum) dot = -umin*sum; 562 #else 563 else if (dot < umin*sum && dot >= 0.0) dot = umin*sum; 564 else if (dot < 0.0 && dot > -umin*sum) dot = -umin*sum; 565 #endif 566 *h = ctx->error_rel*dot/(norm*norm); 567 PetscFunctionReturn(0); 568 } 569