Actual source code: tsfd.c

  1: #define PETSCTS_DLL

 3:  #include private/tsimpl.h

  7: /*@C
  8:     TSDefaultComputeJacobianColor - Computes the Jacobian using
  9:     finite differences and coloring to exploit matrix sparsity.  
 10:   
 11:     Collective on TS, Vec and Mat

 13:     Input Parameters:
 14: +   ts - nonlinear solver object
 15: .   t - current time
 16: .   x1 - location at which to evaluate Jacobian
 17: -   ctx - coloring context, where ctx must have type MatFDColoring, 
 18:           as created via MatFDColoringCreate()

 20:     Output Parameters:
 21: +   J - Jacobian matrix (not altered in this routine)
 22: .   B - newly computed Jacobian matrix to use with preconditioner (generally the same as J)
 23: -   flag - flag indicating whether the matrix sparsity structure has changed

 25:    Level: intermediate

 27: .keywords: TS, finite differences, Jacobian, coloring, sparse

 29: .seealso: TSSetJacobian(), MatFDColoringCreate(), MatFDColoringSetFunction()
 30: @*/
 31: PetscErrorCode  TSDefaultComputeJacobianColor(TS ts,PetscReal t,Vec x1,Mat *J,Mat *B,MatStructure *flag,void *ctx)
 32: {
 33:   MatFDColoring  color = (MatFDColoring) ctx;

 37:   MatFDColoringApplyTS(*B,color,t,x1,flag,ts);
 38: 
 39:   if (*J != *B) {
 40:     MatAssemblyBegin(*J,MAT_FINAL_ASSEMBLY);
 41:     MatAssemblyEnd(*J,MAT_FINAL_ASSEMBLY);
 42:   }
 43:   return(0);
 44: }

 48: /*@C
 49:    TSDefaultComputeJacobian - Computes the Jacobian using finite differences.

 51:    Input Parameters:
 52: +  ts - TS context
 53: .  xx1 - compute Jacobian at this point
 54: -  ctx - application's function context, as set with SNESSetFunction()

 56:    Output Parameters:
 57: +  J - Jacobian
 58: .  B - newly computed Jacobian matrix to use with preconditioner (generally the same as J)
 59: -  flag - matrix flag

 61:    Notes:
 62:    This routine is slow and expensive, and is not optimized.

 64:    Sparse approximations using colorings are also available and
 65:    would be a much better alternative!

 67:    Level: intermediate

 69: .seealso: TSDefaultComputeJacobianColor()
 70: @*/
 71: PetscErrorCode TSDefaultComputeJacobian(TS ts,PetscReal t,Vec xx1,Mat *J,Mat *B,MatStructure *flag,void *ctx)
 72: {
 73:   Vec            f1,f2,xx2;
 75:   PetscInt       i,N,start,end,j;
 76:   PetscScalar    dx,*y,*xx,wscale;
 77:   PetscReal      amax,epsilon = PETSC_SQRT_MACHINE_EPSILON;
 78:   PetscReal      dx_min = 1.e-16,dx_par = 1.e-1;
 79:   MPI_Comm       comm;
 80:   PetscTruth     assembled;
 81:   PetscMPIInt    size;
 82:   const PetscInt *ranges;
 83:   PetscMPIInt    root;

 86:   VecDuplicate(xx1,&f1);
 87:   VecDuplicate(xx1,&f2);
 88:   VecDuplicate(xx1,&xx2);

 90:   PetscObjectGetComm((PetscObject)xx1,&comm);
 91:   MPI_Comm_size(comm,&size);
 92:   MatAssembled(*B,&assembled);
 93:   if (assembled) {
 94:     MatZeroEntries(*B);
 95:   }

 97:   VecGetSize(xx1,&N);
 98:   VecGetOwnershipRange(xx1,&start,&end);
 99:   TSComputeRHSFunction(ts,ts->ptime,xx1,f1);

101:   /* Compute Jacobian approximation, 1 column at a time.
102:       xx1 = current iterate, f1 = F(xx1)
103:       xx2 = perturbed iterate, f2 = F(xx2)
104:    */
105:   for (i=0; i<N; i++) {
106:     VecCopy(xx1,xx2);
107:     if (i>= start && i<end) {
108:        VecGetArray(xx1,&xx);
109:       dx   = xx[i-start];
110:        VecRestoreArray(xx1,&xx);
111: #if !defined(PETSC_USE_COMPLEX)
112:       if (dx < dx_min && dx >= 0.0) dx = dx_par;
113:       else if (dx < 0.0 && dx > -dx_min) dx = -dx_par;
114: #else
115:       if (PetscAbsScalar(dx) < dx_min && PetscRealPart(dx) >= 0.0) dx = dx_par;
116:       else if (PetscRealPart(dx) < 0.0 && PetscAbsScalar(dx) < dx_min) dx = -dx_par;
117: #endif
118:       dx *= epsilon;
119:       wscale = 1.0/dx;
120:        VecSetValues(xx2,1,&i,&dx,ADD_VALUES);
121:     } else {
122:       wscale = 0.0;
123:     }
124:     TSComputeRHSFunction(ts,t,xx2,f2);
125:     VecAXPY(f2,-1.0,f1);
126:     /* Communicate scale=1/dx_i to all processors */
127:     VecGetOwnershipRanges(xx1,&ranges);
128:     root = size;
129:     for (j=size-1; j>-1; j--){
130:       root--;
131:       if (i>=ranges[j]) break;
132:     }
133:     MPI_Bcast(&wscale,1,MPIU_SCALAR,root,comm);

135:     VecScale(f2,wscale);
136:     VecNorm(f2,NORM_INFINITY,&amax); amax *= 1.e-14;
137:     VecGetArray(f2,&y);
138:     for (j=start; j<end; j++) {
139:       if (PetscAbsScalar(y[j-start]) > amax || j == i) {
140:         MatSetValues(*B,1,&j,1,&i,y+j-start,INSERT_VALUES);
141:       }
142:     }
143:     VecRestoreArray(f2,&y);
144:   }
145:   MatAssemblyBegin(*B,MAT_FINAL_ASSEMBLY);
146:   MatAssemblyEnd(*B,MAT_FINAL_ASSEMBLY);
147:   if (*B != *J) {
148:     MatAssemblyBegin(*J,MAT_FINAL_ASSEMBLY);
149:     MatAssemblyEnd(*J,MAT_FINAL_ASSEMBLY);
150:   }
151:   *flag =  DIFFERENT_NONZERO_PATTERN;

153:   VecDestroy(f1);
154:   VecDestroy(f2);
155:   VecDestroy(xx2);
156:   return(0);
157: }