Actual source code: tsfd.c

  1: /*$Id: tsfd.c,v 1.25 2001/08/21 21:04:19 bsmith Exp $*/

 3:  #include src/mat/matimpl.h
 4:  #include src/ts/tsimpl.h

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

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

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

 26:    Options Database Keys:
 27: $  -mat_fd_coloring_freq <freq> 

 29:    Level: intermediate

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

 33: .seealso: TSSetJacobian(), MatFDColoringCreate(), MatFDColoringSetFunction()
 34: @*/
 35: int TSDefaultComputeJacobianColor(TS ts,PetscReal t,Vec x1,Mat *J,Mat *B,MatStructure *flag,void *ctx)
 36: {
 37:   MatFDColoring color = (MatFDColoring) ctx;
 38:   SNES          snes;
 39:   int           ierr,freq,it;

 42:   /*
 43:        If we are not using SNES we have no way to know the current iteration.
 44:   */
 45:   TSGetSNES(ts,&snes);
 46:   if (snes) {
 47:     MatFDColoringGetFrequency(color,&freq);
 48:     SNESGetIterationNumber(snes,&it);

 50:     if ((freq > 1) && ((it % freq) != 1)) {
 51:       PetscLogInfo(color,"TSDefaultComputeJacobianColor:Skipping Jacobian, it %d, freq %dn",it,freq);
 52:       *flag = SAME_PRECONDITIONER;
 53:       return(0);
 54:     } else {
 55:       PetscLogInfo(color,"TSDefaultComputeJacobianColor:Computing Jacobian, it %d, freq %dn",it,freq);
 56:       *flag = SAME_NONZERO_PATTERN;
 57:     }
 58:   }
 59:   MatFDColoringApplyTS(*B,color,t,x1,flag,ts);
 60:   return(0);
 61: }

 63: #undef __FUNCT__  
 65: /*
 66:    TSDefaultComputeJacobian - Computes the Jacobian using finite differences.

 68:    Input Parameters:
 69: +  ts - TS context
 70: .  xx1 - compute Jacobian at this point
 71: -  ctx - application's function context, as set with SNESSetFunction()

 73:    Output Parameters:
 74: +  J - Jacobian
 75: .  B - preconditioner, same as Jacobian
 76: -  flag - matrix flag

 78:    Notes:
 79:    This routine is slow and expensive, and is not optimized.

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

 84:    Level: intermediate

 86: .seealso: TSDefaultComputeJacobianColor()
 87: */
 88: int TSDefaultComputeJacobian(TS ts,PetscReal t,Vec xx1,Mat *J,Mat *B,MatStructure *flag,void *ctx)
 89: {
 90:   Vec         jj1,jj2,xx2;
 91:   int         i,ierr,N,start,end,j;
 92:   PetscScalar dx,mone = -1.0,*y,scale,*xx,wscale;
 93:   PetscReal   amax,epsilon = PETSC_SQRT_MACHINE_EPSILON;
 94:   PetscReal   dx_min = 1.e-16,dx_par = 1.e-1;
 95:   MPI_Comm    comm;

 98:   VecDuplicate(xx1,&jj1);
 99:   VecDuplicate(xx1,&jj2);
100:   VecDuplicate(xx1,&xx2);

102:   PetscObjectGetComm((PetscObject)xx1,&comm);
103:   MatZeroEntries(*J);

105:   VecGetSize(xx1,&N);
106:   VecGetOwnershipRange(xx1,&start,&end);
107:   TSComputeRHSFunction(ts,ts->ptime,xx1,jj1);

109:   /* Compute Jacobian approximation, 1 column at a time.
110:       xx1 = current iterate, jj1 = F(xx1)
111:       xx2 = perturbed iterate, jj2 = F(xx2)
112:    */
113:   for (i=0; i<N; i++) {
114:     VecCopy(xx1,xx2);
115:     if (i>= start && i<end) {
116:        VecGetArray(xx1,&xx);
117:       dx   = xx[i-start];
118:        VecRestoreArray(xx1,&xx);
119: #if !defined(PETSC_USE_COMPLEX)
120:       if (dx < dx_min && dx >= 0.0) dx = dx_par;
121:       else if (dx < 0.0 && dx > -dx_min) dx = -dx_par;
122: #else
123:       if (PetscAbsScalar(dx) < dx_min && PetscRealPart(dx) >= 0.0) dx = dx_par;
124:       else if (PetscRealPart(dx) < 0.0 && PetscAbsScalar(dx) < dx_min) dx = -dx_par;
125: #endif
126:       dx *= epsilon;
127:       wscale = 1.0/dx;
128:        VecSetValues(xx2,1,&i,&dx,ADD_VALUES);
129:     } else {
130:       wscale = 0.0;
131:     }
132:     TSComputeRHSFunction(ts,t,xx2,jj2);
133:     VecAXPY(&mone,jj1,jj2);
134:     /* Communicate scale to all processors */
135:     MPI_Allreduce(&wscale,&scale,1,MPIU_SCALAR,PetscSum_Op,comm);
136:     VecScale(&scale,jj2);
137:     VecNorm(jj2,NORM_INFINITY,&amax); amax *= 1.e-14;
138:     VecGetArray(jj2,&y);
139:     for (j=start; j<end; j++) {
140:       if (PetscAbsScalar(y[j-start]) > amax) {
141:         MatSetValues(*J,1,&j,1,&i,y+j-start,INSERT_VALUES);
142:       }
143:     }
144:     VecRestoreArray(jj2,&y);
145:   }
146:   MatAssemblyBegin(*J,MAT_FINAL_ASSEMBLY);
147:   MatAssemblyEnd(*J,MAT_FINAL_ASSEMBLY);
148:   *flag =  DIFFERENT_NONZERO_PATTERN;

150:   VecDestroy(jj1);
151:   VecDestroy(jj2);
152:   VecDestroy(xx2);

154:   return(0);
155: }