1 #define PETSCMAT_DLL 2 3 /*MC 4 MATMFFD_WP - Implements an alternative approach for computing the differencing parameter 5 h used with the finite difference based matrix-free Jacobian. This code 6 implements the strategy of M. Pernice and H. Walker: 7 8 h = error_rel * sqrt(1 + ||U||) / ||a|| 9 10 Notes: 11 1) || U || does not change between linear iterations so is reused 12 2) In GMRES || a || == 1 and so does not need to ever be computed except at restart 13 when it is recomputed. 14 15 Reference: M. Pernice and H. F. Walker, "NITSOL: A Newton Iterative 16 Solver for Nonlinear Systems", SIAM J. Sci. Stat. Comput.", 1998, 17 vol 19, pp. 302--318. 18 19 Options Database Keys: 20 . -mat_mffd_compute_normu -Compute the norm of u everytime see MatMFFDWPSetComputeNormU() 21 22 23 Level: intermediate 24 25 Notes: 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 "include/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 PetscTruth computenorma,computenormU; 47 } MatMFFD_WP; 48 49 #undef __FUNCT__ 50 #define __FUNCT__ "MatMFFDCompute_WP" 51 /* 52 MatMFFDCompute_WP - Standard PETSc code for 53 computing h with matrix-free finite differences. 54 55 Input Parameters: 56 + ctx - the matrix free context 57 . U - the location at which you want the Jacobian 58 - a - the direction you want the derivative 59 60 Output Parameter: 61 . h - the scale computed 62 63 */ 64 static PetscErrorCode MatMFFDCompute_WP(MatMFFD ctx,Vec U,Vec a,PetscScalar *h,PetscTruth *zeroa) 65 { 66 MatMFFD_WP *hctx = (MatMFFD_WP*)ctx->hctx; 67 PetscReal normU,norma = 1.0; 68 PetscErrorCode ierr; 69 70 PetscFunctionBegin; 71 if (!(ctx->count % ctx->recomputeperiod)) { 72 if (hctx->computenorma && (hctx->computenormU || !ctx->ncurrenth)) { 73 ierr = VecNormBegin(U,NORM_2,&normU);CHKERRQ(ierr); 74 ierr = VecNormBegin(a,NORM_2,&norma);CHKERRQ(ierr); 75 ierr = VecNormEnd(U,NORM_2,&normU);CHKERRQ(ierr); 76 ierr = VecNormEnd(a,NORM_2,&norma);CHKERRQ(ierr); 77 hctx->normUfact = sqrt(1.0+normU); 78 } else if (hctx->computenormU || !ctx->ncurrenth) { 79 ierr = VecNorm(U,NORM_2,&normU);CHKERRQ(ierr); 80 hctx->normUfact = sqrt(1.0+normU); 81 } else if (hctx->computenorma) { 82 ierr = VecNorm(a,NORM_2,&norma);CHKERRQ(ierr); 83 } 84 if (norma == 0.0) { 85 *zeroa = PETSC_TRUE; 86 PetscFunctionReturn(0); 87 } 88 *zeroa = PETSC_FALSE; 89 *h = ctx->error_rel*hctx->normUfact/norma; 90 } else { 91 *h = ctx->currenth; 92 } 93 PetscFunctionReturn(0); 94 } 95 96 #undef __FUNCT__ 97 #define __FUNCT__ "MatMFFDView_WP" 98 /* 99 MatMFFDView_WP - Prints information about this particular 100 method for computing h. Note that this does not print the general 101 information about the matrix free, that is printed by the calling 102 routine. 103 104 Input Parameters: 105 + ctx - the matrix free context 106 - viewer - the PETSc viewer 107 108 */ 109 static PetscErrorCode MatMFFDView_WP(MatMFFD ctx,PetscViewer viewer) 110 { 111 MatMFFD_WP *hctx = (MatMFFD_WP *)ctx->hctx; 112 PetscErrorCode ierr; 113 PetscTruth iascii; 114 115 PetscFunctionBegin; 116 ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr); 117 if (iascii) { 118 if (hctx->computenorma){ierr = PetscViewerASCIIPrintf(viewer," Computes normA\n");CHKERRQ(ierr);} 119 else {ierr = PetscViewerASCIIPrintf(viewer," Does not compute normA\n");CHKERRQ(ierr);} 120 if (hctx->computenormU){ierr = PetscViewerASCIIPrintf(viewer," Computes normU\n");CHKERRQ(ierr);} 121 else {ierr = PetscViewerASCIIPrintf(viewer," Does not compute normU\n");CHKERRQ(ierr);} 122 } else { 123 SETERRQ1(PETSC_ERR_SUP,"Viewer type %s not supported for SNES matrix-free WP",((PetscObject)viewer)->type_name); 124 } 125 PetscFunctionReturn(0); 126 } 127 128 #undef __FUNCT__ 129 #define __FUNCT__ "MatMFFDSetFromOptions_WP" 130 /* 131 MatMFFDSetFromOptions_WP - Looks in the options database for 132 any options appropriate for this method 133 134 Input Parameter: 135 . ctx - the matrix free context 136 137 */ 138 static PetscErrorCode MatMFFDSetFromOptions_WP(MatMFFD ctx) 139 { 140 PetscErrorCode ierr; 141 MatMFFD_WP *hctx = (MatMFFD_WP*)ctx->hctx; 142 143 PetscFunctionBegin; 144 ierr = PetscOptionsHead("Walker-Pernice options");CHKERRQ(ierr); 145 ierr = PetscOptionsTruth("-mat_mffd_compute_normu","Compute the norm of u","MatMFFDWPSetComputeNormU", 146 hctx->computenorma,&hctx->computenorma,0);CHKERRQ(ierr); 147 ierr = PetscOptionsTail();CHKERRQ(ierr); 148 PetscFunctionReturn(0); 149 } 150 151 #undef __FUNCT__ 152 #define __FUNCT__ "MatMFFDDestroy_WP" 153 /* 154 MatMFFDDestroy_WP - Frees the space allocated by 155 MatMFFDCreate_WP(). 156 157 Input Parameter: 158 . ctx - the matrix free context 159 160 Notes: does not free the ctx, that is handled by the calling routine 161 162 */ 163 static PetscErrorCode MatMFFDDestroy_WP(MatMFFD ctx) 164 { 165 PetscErrorCode ierr; 166 PetscFunctionBegin; 167 ierr = PetscFree(ctx->hctx);CHKERRQ(ierr); 168 PetscFunctionReturn(0); 169 } 170 171 EXTERN_C_BEGIN 172 #undef __FUNCT__ 173 #define __FUNCT__ "MatMFFDWPSetComputeNormU_P" 174 PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDWPSetComputeNormU_P(Mat mat,PetscTruth flag) 175 { 176 MatMFFD ctx = (MatMFFD)mat->data; 177 MatMFFD_WP *hctx = (MatMFFD_WP*)ctx->hctx; 178 179 PetscFunctionBegin; 180 hctx->computenormU = flag; 181 PetscFunctionReturn(0); 182 } 183 EXTERN_C_END 184 185 #undef __FUNCT__ 186 #define __FUNCT__ "MatMFFDWPSetComputeNormU" 187 /*@ 188 MatMFFDWPSetComputeNormU - Sets whether it computes the ||U|| used by the WP 189 PETSc routine for computing h. With any Krylov solver this need only 190 be computed during the first iteration and kept for later. 191 192 Input Parameters: 193 + A - the matrix created with MatCreateSNESMF() 194 - flag - PETSC_TRUE causes it to compute ||U||, PETSC_FALSE uses the previous value 195 196 Options Database Key: 197 . -mat_mffd_compute_normu <true,false> - true by default, false can save calculations but you 198 must be sure that ||U|| has not changed in the mean time. 199 200 Level: advanced 201 202 Notes: 203 See the manual page for MATMFFD_WP for a complete description of the 204 algorithm used to compute h. 205 206 .seealso: MatMFFDSetFunctionError(), MatCreateSNESMF() 207 208 @*/ 209 PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDWPSetComputeNormU(Mat A,PetscTruth flag) 210 { 211 PetscErrorCode ierr,(*f)(Mat,PetscTruth); 212 213 PetscFunctionBegin; 214 PetscValidHeaderSpecific(A,MAT_COOKIE,1); 215 ierr = PetscObjectQueryFunction((PetscObject)A,"MatMFFDWPSetComputeNormU_C",(void (**)(void))&f);CHKERRQ(ierr); 216 if (f) { 217 ierr = (*f)(A,flag);CHKERRQ(ierr); 218 } 219 PetscFunctionReturn(0); 220 } 221 222 EXTERN_C_BEGIN 223 #undef __FUNCT__ 224 #define __FUNCT__ "MatMFFDCreate_WP" 225 /* 226 MatMFFDCreate_WP - Standard PETSc code for 227 computing h with matrix-free finite differences. 228 229 Input Parameter: 230 . ctx - the matrix free context created by MatMFFDCreate() 231 232 */ 233 PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDCreate_WP(MatMFFD ctx) 234 { 235 PetscErrorCode ierr; 236 MatMFFD_WP *hctx; 237 238 PetscFunctionBegin; 239 240 /* allocate my own private data structure */ 241 ierr = PetscNewLog(ctx,MatMFFD_WP,&hctx);CHKERRQ(ierr); 242 ctx->hctx = (void*)hctx; 243 hctx->computenormU = PETSC_FALSE; 244 hctx->computenorma = PETSC_TRUE; 245 246 /* set the functions I am providing */ 247 ctx->ops->compute = MatMFFDCompute_WP; 248 ctx->ops->destroy = MatMFFDDestroy_WP; 249 ctx->ops->view = MatMFFDView_WP; 250 ctx->ops->setfromoptions = MatMFFDSetFromOptions_WP; 251 252 ierr = PetscObjectComposeFunctionDynamic((PetscObject)ctx->mat,"MatMFFDWPSetComputeNormU_C", 253 "MatMFFDWPSetComputeNormU_P", 254 MatMFFDWPSetComputeNormU_P);CHKERRQ(ierr); 255 256 PetscFunctionReturn(0); 257 } 258 EXTERN_C_END 259 260 261 262