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: /*@C
7: TSDefaultComputeJacobianColor - Computes the Jacobian using
8: finite differences and coloring to exploit matrix sparsity.
9:
10: Collective on TS, Vec and Mat
12: Input Parameters:
13: + ts - nonlinear solver object
14: . t - current time
15: . x1 - location at which to evaluate Jacobian
16: - ctx - coloring context, where ctx must have type MatFDColoring,
17: as created via MatFDColoringCreate()
19: Output Parameters:
20: + J - Jacobian matrix (not altered in this routine)
21: . B - newly computed Jacobian matrix to use with preconditioner (generally the same as J)
22: - flag - flag indicating whether the matrix sparsity structure has changed
24: Options Database Keys:
25: $ -mat_fd_coloring_freq <freq>
27: Level: intermediate
29: .keywords: TS, finite differences, Jacobian, coloring, sparse
31: .seealso: TSSetJacobian(), MatFDColoringCreate(), MatFDColoringSetFunction()
32: @*/
33: int TSDefaultComputeJacobianColor(TS ts,PetscReal t,Vec x1,Mat *J,Mat *B,MatStructure *flag,void *ctx)
34: {
35: MatFDColoring color = (MatFDColoring) ctx;
36: SNES snes;
37: int ierr,freq,it;
40: /*
41: If we are not using SNES we have no way to know the current iteration.
42: */
43: TSGetSNES(ts,&snes);
44: if (snes) {
45: MatFDColoringGetFrequency(color,&freq);
46: SNESGetIterationNumber(snes,&it);
48: if ((freq > 1) && ((it % freq) != 1)) {
49: PetscLogInfo(color,"TSDefaultComputeJacobianColor:Skipping Jacobian, it %d, freq %dn",it,freq);
50: *flag = SAME_PRECONDITIONER;
51: return(0);
52: } else {
53: PetscLogInfo(color,"TSDefaultComputeJacobianColor:Computing Jacobian, it %d, freq %dn",it,freq);
54: *flag = SAME_NONZERO_PATTERN;
55: }
56: }
57: MatFDColoringApplyTS(*B,color,t,x1,flag,ts);
58: return(0);
59: }
61: /*
62: TSDefaultComputeJacobian - Computes the Jacobian using finite differences.
64: Input Parameters:
65: + ts - TS context
66: . xx1 - compute Jacobian at this point
67: - ctx - application's function context, as set with SNESSetFunction()
69: Output Parameters:
70: + J - Jacobian
71: . B - preconditioner, same as Jacobian
72: - flag - matrix flag
74: Notes:
75: This routine is slow and expensive, and is not optimized.
77: Sparse approximations using colorings are also available and
78: would be a much better alternative!
80: Level: intermediate
82: .seealso: TSDefaultComputeJacobianColor()
83: */
84: int TSDefaultComputeJacobian(TS ts,PetscReal t,Vec xx1,Mat *J,Mat *B,MatStructure *flag,void *ctx)
85: {
86: Vec jj1,jj2,xx2;
87: int i,ierr,N,start,end,j;
88: PetscScalar dx,mone = -1.0,*y,scale,*xx,wscale;
89: PetscReal amax,epsilon = PETSC_SQRT_MACHINE_EPSILON;
90: PetscReal dx_min = 1.e-16,dx_par = 1.e-1;
91: MPI_Comm comm;
94: VecDuplicate(xx1,&jj1);
95: VecDuplicate(xx1,&jj2);
96: VecDuplicate(xx1,&xx2);
98: PetscObjectGetComm((PetscObject)xx1,&comm);
99: MatZeroEntries(*J);
101: VecGetSize(xx1,&N);
102: VecGetOwnershipRange(xx1,&start,&end);
103: TSComputeRHSFunction(ts,ts->ptime,xx1,jj1);
105: /* Compute Jacobian approximation, 1 column at a time.
106: xx1 = current iterate, jj1 = F(xx1)
107: xx2 = perturbed iterate, jj2 = F(xx2)
108: */
109: for (i=0; i<N; i++) {
110: VecCopy(xx1,xx2);
111: if (i>= start && i<end) {
112: VecGetArray(xx1,&xx);
113: dx = xx[i-start];
114: VecRestoreArray(xx1,&xx);
115: #if !defined(PETSC_USE_COMPLEX)
116: if (dx < dx_min && dx >= 0.0) dx = dx_par;
117: else if (dx < 0.0 && dx > -dx_min) dx = -dx_par;
118: #else
119: if (PetscAbsScalar(dx) < dx_min && PetscRealPart(dx) >= 0.0) dx = dx_par;
120: else if (PetscRealPart(dx) < 0.0 && PetscAbsScalar(dx) < dx_min) dx = -dx_par;
121: #endif
122: dx *= epsilon;
123: wscale = 1.0/dx;
124: VecSetValues(xx2,1,&i,&dx,ADD_VALUES);
125: } else {
126: wscale = 0.0;
127: }
128: TSComputeRHSFunction(ts,t,xx2,jj2);
129: VecAXPY(&mone,jj1,jj2);
130: /* Communicate scale to all processors */
131: MPI_Allreduce(&wscale,&scale,1,MPIU_SCALAR,PetscSum_Op,comm);
132: VecScale(&scale,jj2);
133: VecNorm(jj2,NORM_INFINITY,&amax); amax *= 1.e-14;
134: VecGetArray(jj2,&y);
135: for (j=start; j<end; j++) {
136: if (PetscAbsScalar(y[j-start]) > amax) {
137: MatSetValues(*J,1,&j,1,&i,y+j-start,INSERT_VALUES);
138: }
139: }
140: VecRestoreArray(jj2,&y);
141: }
142: MatAssemblyBegin(*J,MAT_FINAL_ASSEMBLY);
143: MatAssemblyEnd(*J,MAT_FINAL_ASSEMBLY);
144: *flag = DIFFERENT_NONZERO_PATTERN;
146: VecDestroy(jj1);
147: VecDestroy(jj2);
148: VecDestroy(xx2);
150: return(0);
151: }