1 2 /*MC 3 MATMFFD_WP - Implements an alternative approach for computing the differencing parameter 4 h used with the finite difference based matrix-free Jacobian. This code 5 implements the strategy of M. Pernice and H. Walker: 6 7 h = error_rel * sqrt(1 + ||U||) / ||a|| 8 9 Notes: 10 1) || U || does not change between linear iterations so is reused 11 2) In GMRES || a || == 1 and so does not need to ever be computed except at restart 12 when it is recomputed. 13 14 Reference: M. Pernice and H. F. Walker, "NITSOL: A Newton Iterative 15 Solver for Nonlinear Systems", SIAM J. Sci. Stat. Comput.", 1998, 16 vol 19, pp. 302--318. 17 18 Options Database Keys: 19 . -mat_mffd_compute_normu -Compute the norm of u everytime see MatMFFDWPSetComputeNormU() 20 21 22 Level: intermediate 23 24 Notes: 25 Requires no global collectives when used with GMRES 26 27 Formula used: 28 F'(u)*a = [F(u+h*a) - F(u)]/h where 29 30 .seealso: MATMFFD, MatCreateMFFD(), MatCreateSNESMF(), MATMFFD_DS 31 32 M*/ 33 34 /* 35 This include file defines the data structure MatMFFD that 36 includes information about the computation of h. It is shared by 37 all implementations that people provide. 38 39 See snesmfjdef.c for a full set of comments on the routines below. 40 */ 41 #include <petsc/private/matimpl.h> 42 #include <../src/mat/impls/mffd/mffdimpl.h> /*I "petscmat.h" I*/ 43 44 typedef struct { 45 PetscReal normUfact; /* previous sqrt(1.0 + || U ||) */ 46 PetscBool computenormU; 47 } MatMFFD_WP; 48 49 /* 50 MatMFFDCompute_WP - Standard PETSc code for 51 computing h with matrix-free finite differences. 52 53 Input Parameters: 54 + ctx - the matrix free context 55 . U - the location at which you want the Jacobian 56 - a - the direction you want the derivative 57 58 Output Parameter: 59 . h - the scale computed 60 61 */ 62 static PetscErrorCode MatMFFDCompute_WP(MatMFFD ctx,Vec U,Vec a,PetscScalar *h,PetscBool *zeroa) 63 { 64 MatMFFD_WP *hctx = (MatMFFD_WP*)ctx->hctx; 65 PetscReal normU,norma; 66 PetscErrorCode ierr; 67 68 PetscFunctionBegin; 69 if (!(ctx->count % ctx->recomputeperiod)) { 70 if (hctx->computenormU || !ctx->ncurrenth) { 71 ierr = VecNorm(U,NORM_2,&normU);CHKERRQ(ierr); 72 hctx->normUfact = PetscSqrtReal(1.0+normU); 73 } 74 ierr = VecNorm(a,NORM_2,&norma);CHKERRQ(ierr); 75 if (norma == 0.0) { 76 *zeroa = PETSC_TRUE; 77 PetscFunctionReturn(0); 78 } 79 *zeroa = PETSC_FALSE; 80 *h = ctx->error_rel*hctx->normUfact/norma; 81 } else { 82 *h = ctx->currenth; 83 } 84 PetscFunctionReturn(0); 85 } 86 87 /* 88 MatMFFDView_WP - Prints information about this particular 89 method for computing h. Note that this does not print the general 90 information about the matrix free, that is printed by the calling 91 routine. 92 93 Input Parameters: 94 + ctx - the matrix free context 95 - viewer - the PETSc viewer 96 97 */ 98 static PetscErrorCode MatMFFDView_WP(MatMFFD ctx,PetscViewer viewer) 99 { 100 MatMFFD_WP *hctx = (MatMFFD_WP*)ctx->hctx; 101 PetscErrorCode ierr; 102 PetscBool iascii; 103 104 PetscFunctionBegin; 105 ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 106 if (iascii) { 107 if (hctx->computenormU) { 108 ierr = PetscViewerASCIIPrintf(viewer," Computes normU\n");CHKERRQ(ierr); 109 } else { 110 ierr = PetscViewerASCIIPrintf(viewer," Does not compute normU\n");CHKERRQ(ierr); 111 } 112 } 113 PetscFunctionReturn(0); 114 } 115 116 /* 117 MatMFFDSetFromOptions_WP - Looks in the options database for 118 any options appropriate for this method 119 120 Input Parameter: 121 . ctx - the matrix free context 122 123 */ 124 static PetscErrorCode MatMFFDSetFromOptions_WP(PetscOptionItems *PetscOptionsObject,MatMFFD ctx) 125 { 126 PetscErrorCode ierr; 127 MatMFFD_WP *hctx = (MatMFFD_WP*)ctx->hctx; 128 129 PetscFunctionBegin; 130 ierr = PetscOptionsHead(PetscOptionsObject,"Walker-Pernice options");CHKERRQ(ierr); 131 ierr = PetscOptionsBool("-mat_mffd_compute_normu","Compute the norm of u","MatMFFDWPSetComputeNormU", hctx->computenormU,&hctx->computenormU,NULL);CHKERRQ(ierr); 132 ierr = PetscOptionsTail();CHKERRQ(ierr); 133 PetscFunctionReturn(0); 134 } 135 136 /* 137 MatMFFDDestroy_WP - Frees the space allocated by 138 MatCreateMFFD_WP(). 139 140 Input Parameter: 141 . ctx - the matrix free context 142 143 Notes: 144 does not free the ctx, that is handled by the calling routine 145 146 */ 147 static PetscErrorCode MatMFFDDestroy_WP(MatMFFD ctx) 148 { 149 PetscErrorCode ierr; 150 151 PetscFunctionBegin; 152 ierr = PetscFree(ctx->hctx);CHKERRQ(ierr); 153 PetscFunctionReturn(0); 154 } 155 156 PetscErrorCode MatMFFDWPSetComputeNormU_P(Mat mat,PetscBool flag) 157 { 158 MatMFFD ctx = (MatMFFD)mat->data; 159 MatMFFD_WP *hctx = (MatMFFD_WP*)ctx->hctx; 160 161 PetscFunctionBegin; 162 hctx->computenormU = flag; 163 PetscFunctionReturn(0); 164 } 165 166 /*@ 167 MatMFFDWPSetComputeNormU - Sets whether it computes the ||U|| used by the WP 168 PETSc routine for computing h. With any Krylov solver this need only 169 be computed during the first iteration and kept for later. 170 171 Input Parameters: 172 + A - the matrix created with MatCreateSNESMF() 173 - flag - PETSC_TRUE causes it to compute ||U||, PETSC_FALSE uses the previous value 174 175 Options Database Key: 176 . -mat_mffd_compute_normu <true,false> - true by default, false can save calculations but you 177 must be sure that ||U|| has not changed in the mean time. 178 179 Level: advanced 180 181 Notes: 182 See the manual page for MATMFFD_WP for a complete description of the 183 algorithm used to compute h. 184 185 .seealso: MatMFFDSetFunctionError(), MatCreateSNESMF() 186 187 @*/ 188 PetscErrorCode MatMFFDWPSetComputeNormU(Mat A,PetscBool flag) 189 { 190 PetscErrorCode ierr; 191 192 PetscFunctionBegin; 193 PetscValidHeaderSpecific(A,MAT_CLASSID,1); 194 ierr = PetscTryMethod(A,"MatMFFDWPSetComputeNormU_C",(Mat,PetscBool),(A,flag));CHKERRQ(ierr); 195 PetscFunctionReturn(0); 196 } 197 198 /* 199 MatCreateMFFD_WP - Standard PETSc code for 200 computing h with matrix-free finite differences. 201 202 Input Parameter: 203 . ctx - the matrix free context created by MatCreateMFFD() 204 205 */ 206 PETSC_EXTERN PetscErrorCode MatCreateMFFD_WP(MatMFFD ctx) 207 { 208 PetscErrorCode ierr; 209 MatMFFD_WP *hctx; 210 211 PetscFunctionBegin; 212 /* allocate my own private data structure */ 213 ierr = PetscNewLog(ctx,&hctx);CHKERRQ(ierr); 214 ctx->hctx = (void*)hctx; 215 hctx->computenormU = PETSC_FALSE; 216 217 /* set the functions I am providing */ 218 ctx->ops->compute = MatMFFDCompute_WP; 219 ctx->ops->destroy = MatMFFDDestroy_WP; 220 ctx->ops->view = MatMFFDView_WP; 221 ctx->ops->setfromoptions = MatMFFDSetFromOptions_WP; 222 223 ierr = PetscObjectComposeFunction((PetscObject)ctx->mat,"MatMFFDWPSetComputeNormU_C",MatMFFDWPSetComputeNormU_P);CHKERRQ(ierr); 224 PetscFunctionReturn(0); 225 } 226 227 228 229