Subversion Repositories slepc-dev

Rev

Show entire file | Ignore whitespace | Details | Blame | Last modification | View Log | RSS feed

Rev 1229 Rev 1345
Line 204... Line 204...
  ierr = PetscFree(rows);CHKERRQ(ierr);
  ierr = PetscFree(rows);CHKERRQ(ierr);
  ierr = VecDestroy(in);CHKERRQ(ierr);
  ierr = VecDestroy(in);CHKERRQ(ierr);
  ierr = VecDestroy(out);CHKERRQ(ierr);
  ierr = VecDestroy(out);CHKERRQ(ierr);
  ierr = MatAssemblyBegin(*mat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
  ierr = MatAssemblyBegin(*mat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
  ierr = MatAssemblyEnd(*mat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
  ierr = MatAssemblyEnd(*mat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
  PetscFunctionReturn(0);
 
}
 
 
 
#undef __FUNCT__  
 
#define __FUNCT__ "STNorm"
 
/*@
 
   STNorm - Computes de norm of a vector as the square root of the inner
 
   product (x,x) as defined by STInnerProduct().
 
 
 
   Collective on ST and Vec
 
 
 
   Input Parameters:
 
+  st - the spectral transformation context
 
-  x  - input vector
 
 
 
   Output Parameter:
 
.  norm - the computed norm
 
 
 
   Notes:
 
   This function will usually compute the 2-norm of a vector, ||x||_2. But
 
   this behaviour may be different if using a non-standard inner product changed
 
   via STSetBilinearForm(). For example, if using the B-inner product for
 
   positive definite B, (x,y)_B=y^H Bx, then the computed norm is ||x||_B =
 
   sqrt( x^H Bx ).
 
 
 
   Level: developer
 
 
 
.seealso: STInnerProduct()
 
@*/
 
PetscErrorCode STNorm(ST st,Vec x,PetscReal *norm)
 
{
 
  PetscErrorCode ierr;
 
  PetscScalar    p;
 
 
 
  PetscFunctionBegin;
 
  PetscValidHeaderSpecific(st,ST_COOKIE,1);
 
  PetscValidHeaderSpecific(x,VEC_COOKIE,2);
 
  PetscValidPointer(norm,3);
 
 
 
  ierr = STInnerProduct(st,x,x,&p);CHKERRQ(ierr);
 
 
 
  if (PetscAbsScalar(p)<PETSC_MACHINE_EPSILON)
 
    PetscInfo(st,"Zero norm, either the vector is zero or a semi-inner product is being used\n");
 
 
 
#if defined(PETSC_USE_COMPLEX)
 
  if (PetscRealPart(p)<0.0 || PetscAbsReal(PetscImaginaryPart(p))>PETSC_MACHINE_EPSILON)
 
     SETERRQ(1,"STNorm: The inner product is not well defined");
 
  *norm = PetscSqrtScalar(PetscRealPart(p));
 
#else
 
  if (p<0.0) SETERRQ(1,"STNorm: The inner product is not well defined");
 
  *norm = PetscSqrtScalar(p);
 
#endif
 
 
 
  PetscFunctionReturn(0);
 
}
 
 
 
#undef __FUNCT__  
 
#define __FUNCT__ "STNormBegin"
 
/*@
 
   STNormBegin - Starts a split phase norm computation.
 
 
 
   Input Parameters:
 
+  st   - the spectral transformation context
 
.  x    - input vector
 
-  norm - where the result will go
 
 
 
   Level: developer
 
 
 
   Notes:
 
   Each call to STNormBegin() should be paired with a call to STNormEnd().
 
 
 
.seealso: STNormEnd(), STNorm(), STInnerProduct(), STMInnerProduct(),
 
          STInnerProductBegin(), STInnerProductEnd()
 
 
 
@*/
 
PetscErrorCode STNormBegin(ST st,Vec x,PetscReal *norm)
 
{
 
  PetscErrorCode ierr;
 
  PetscScalar    p;
 
 
 
  PetscFunctionBegin;
 
  PetscValidHeaderSpecific(st,ST_COOKIE,1);
 
  PetscValidHeaderSpecific(x,VEC_COOKIE,2);
 
  PetscValidPointer(norm,3);
 
 
 
  ierr = STInnerProductBegin(st,x,x,&p);CHKERRQ(ierr);
 
 
 
  PetscFunctionReturn(0);
 
}
 
 
 
#undef __FUNCT__  
 
#define __FUNCT__ "STNormEnd"
 
/*@
 
   STNormEnd - Ends a split phase norm computation.
 
 
 
   Input Parameters:
 
+  st   - the spectral transformation context
 
-  x    - input vector
 
 
 
   Output Parameter:
 
.  norm - the computed norm
 
 
 
   Level: developer
 
 
 
   Notes:
 
   Each call to STNormBegin() should be paired with a call to STNormEnd().
 
 
 
.seealso: STNormBegin(), STNorm(), STInnerProduct(), STMInnerProduct(),
 
          STInnerProductBegin(), STInnerProductEnd()
 
 
 
@*/
 
PetscErrorCode STNormEnd(ST st,Vec x,PetscReal *norm)
 
{
 
  PetscErrorCode ierr;
 
  PetscScalar    p;
 
 
 
  PetscFunctionBegin;
 
  PetscValidHeaderSpecific(st,ST_COOKIE,1);
 
  PetscValidHeaderSpecific(x,VEC_COOKIE,2);
 
  PetscValidPointer(norm,3);
 
 
 
  ierr = STInnerProductEnd(st,x,x,&p);CHKERRQ(ierr);
 
 
 
  if (PetscAbsScalar(p)<PETSC_MACHINE_EPSILON)
 
    PetscInfo(st,"Zero norm, either the vector is zero or a semi-inner product is being used\n");
 
 
 
#if defined(PETSC_USE_COMPLEX)
 
  if (PetscRealPart(p)<0.0 || PetscAbsReal(PetscImaginaryPart(p))>PETSC_MACHINE_EPSILON)
 
     SETERRQ(1,"STNorm: The inner product is not well defined");
 
  *norm = PetscSqrtScalar(PetscRealPart(p));
 
#else
 
  if (p<0.0) SETERRQ(1,"STNorm: The inner product is not well defined");
 
  *norm = PetscSqrtScalar(p);
 
#endif
 
 
 
  PetscFunctionReturn(0);
 
}
 
 
 
#undef __FUNCT__  
 
#define __FUNCT__ "STInnerProduct"
 
/*@
 
   STInnerProduct - Computes the inner product of two vectors.
 
 
 
   Collective on ST and Vec
 
 
 
   Input Parameters:
 
+  st - the spectral transformation context
 
.  x  - input vector
 
-  y  - input vector
 
 
 
   Output Parameter:
 
.  p - result of the inner product
 
 
 
   Notes:
 
   This function will usually compute the standard dot product of vectors
 
   x and y, (x,y)=y^H x. However this behaviour may be different if changed
 
   via STSetBilinearForm(). This allows use of other inner products such as
 
   the indefinite product y^T x for complex symmetric problems or the
 
   B-inner product for positive definite B, (x,y)_B=y^H Bx.
 
 
 
   Level: developer
 
 
 
.seealso: STSetBilinearForm(), STApplyB(), VecDot(), STMInnerProduct()
 
@*/
 
PetscErrorCode STInnerProduct(ST st,Vec x,Vec y,PetscScalar *p)
 
{
 
  PetscErrorCode ierr;
 
 
 
  PetscFunctionBegin;
 
  PetscValidHeaderSpecific(st,ST_COOKIE,1);
 
  PetscValidHeaderSpecific(x,VEC_COOKIE,2);
 
  PetscValidHeaderSpecific(y,VEC_COOKIE,3);
 
  PetscValidScalarPointer(p,4);
 
 
 
  ierr = PetscLogEventBegin(ST_InnerProduct,st,x,0,0);CHKERRQ(ierr);
 
  st->innerproducts++;
 
  switch (st->bilinear_form) {
 
  case STINNER_HERMITIAN:
 
  case STINNER_SYMMETRIC:
 
    ierr = VecCopy(x,st->w);CHKERRQ(ierr);
 
    break;
 
  case STINNER_B_HERMITIAN:
 
  case STINNER_B_SYMMETRIC:
 
    ierr = STApplyB(st,x,st->w);CHKERRQ(ierr);
 
    break;
 
  }
 
  switch (st->bilinear_form) {
 
  case STINNER_HERMITIAN:
 
  case STINNER_B_HERMITIAN:
 
    ierr = VecDot(st->w,y,p);CHKERRQ(ierr);
 
    break;
 
  case STINNER_SYMMETRIC:
 
  case STINNER_B_SYMMETRIC:
 
    ierr = VecTDot(st->w,y,p);CHKERRQ(ierr);
 
    break;
 
  }
 
  ierr = PetscLogEventEnd(ST_InnerProduct,st,x,0,0);CHKERRQ(ierr);
 
  PetscFunctionReturn(0);
 
}
 
 
 
#undef __FUNCT__  
 
#define __FUNCT__ "STInnerProductBegin"
 
/*@
 
   STInnerProductBegin - Starts a split phase inner product computation.
 
 
 
   Input Parameters:
 
+  st - the spectral transformation context
 
.  x  - the first vector
 
.  y  - the second vector
 
-  p  - where the result will go
 
 
 
   Level: developer
 
 
 
   Notes:
 
   Each call to STInnerProductBegin() should be paired with a call to STInnerProductEnd().
 
 
 
.seealso: STInnerProductEnd(), STInnerProduct(), STNorm(), STNormBegin(),
 
          STNormEnd(), STMInnerProduct()
 
 
 
@*/
 
PetscErrorCode STInnerProductBegin(ST st,Vec x,Vec y,PetscScalar *p)
 
{
 
  PetscErrorCode ierr;
 
 
 
  PetscFunctionBegin;
 
  PetscValidHeaderSpecific(st,ST_COOKIE,1);
 
  PetscValidHeaderSpecific(x,VEC_COOKIE,2);
 
  PetscValidHeaderSpecific(y,VEC_COOKIE,3);
 
  PetscValidScalarPointer(p,4);
 
 
 
  ierr = PetscLogEventBegin(ST_InnerProduct,st,x,0,0);CHKERRQ(ierr);
 
  st->innerproducts++;
 
  switch (st->bilinear_form) {
 
  case STINNER_HERMITIAN:
 
  case STINNER_SYMMETRIC:
 
    ierr = VecCopy(x,st->w);CHKERRQ(ierr);
 
    break;
 
  case STINNER_B_HERMITIAN:
 
  case STINNER_B_SYMMETRIC:
 
    ierr = STApplyB(st,x,st->w);CHKERRQ(ierr);
 
    break;
 
  }
 
  switch (st->bilinear_form) {
 
  case STINNER_HERMITIAN:
 
  case STINNER_B_HERMITIAN:
 
    ierr = VecDotBegin(st->w,y,p);CHKERRQ(ierr);
 
    break;
 
  case STINNER_SYMMETRIC:
 
  case STINNER_B_SYMMETRIC:
 
    ierr = VecTDotBegin(st->w,y,p);CHKERRQ(ierr);
 
    break;
 
  }
 
  ierr = PetscLogEventEnd(ST_InnerProduct,st,x,0,0);CHKERRQ(ierr);
 
  PetscFunctionReturn(0);
 
}
 
 
 
#undef __FUNCT__  
 
#define __FUNCT__ "STInnerProductEnd"
 
/*@
 
   STInnerProductEnd - Ends a split phase inner product computation.
 
 
 
   Input Parameters:
 
+  st - the spectral transformation context
 
.  x  - the first vector
 
-  y  - the second vector
 
 
 
   Output Parameter:
 
.  p  - result of the inner product
 
 
 
   Level: developer
 
 
 
   Notes:
 
   Each call to STInnerProductBegin() should be paired with a call to STInnerProductEnd().
 
 
 
.seealso: STInnerProductBegin(), STInnerProduct(), STNorm(), STNormBegin(),
 
          STNormEnd(), STMInnerProduct()
 
 
 
@*/
 
PetscErrorCode STInnerProductEnd(ST st,Vec x,Vec y,PetscScalar *p)
 
{
 
  PetscErrorCode ierr;
 
 
 
  PetscFunctionBegin;
 
  PetscValidHeaderSpecific(st,ST_COOKIE,1);
 
  PetscValidHeaderSpecific(x,VEC_COOKIE,2);
 
  PetscValidHeaderSpecific(y,VEC_COOKIE,3);
 
  PetscValidScalarPointer(p,4);
 
 
 
  ierr = PetscLogEventBegin(ST_InnerProduct,st,x,0,0);CHKERRQ(ierr);
 
  switch (st->bilinear_form) {
 
  case STINNER_HERMITIAN:
 
  case STINNER_B_HERMITIAN:
 
    ierr = VecDotEnd(st->w,y,p);CHKERRQ(ierr);
 
    break;
 
  case STINNER_SYMMETRIC:
 
  case STINNER_B_SYMMETRIC:
 
    ierr = VecTDotEnd(st->w,y,p);CHKERRQ(ierr);
 
    break;
 
  }
 
  ierr = PetscLogEventEnd(ST_InnerProduct,st,x,0,0);CHKERRQ(ierr);
 
  PetscFunctionReturn(0);
 
}
 
 
 
#undef __FUNCT__  
 
#define __FUNCT__ "STMInnerProduct"
 
/*@
 
   STMInnerProduct - Computes the inner products a vector x with a set of
 
   vectors (columns of Y).
 
 
 
   Collective on ST and Vec
 
 
 
   Input Parameters:
 
+  st - the spectral transformation context
 
.  n  - number of vectors in y
 
.  x  - the first input vector
 
-  y  - array of vectors
 
 
 
   Output Parameter:
 
.  p - result of the inner products
 
 
 
   Notes:
 
   This function will usually compute the standard dot product of x and y_i,
 
   (x,y_i)=y_i^H x, for each column of Y. However this behaviour may be different
 
   if changed via STSetBilinearForm(). This allows use of other inner products
 
   such as the indefinite product y_i^T x for complex symmetric problems or the
 
   B-inner product for positive definite B, (x,y_i)_B=y_i^H Bx.
 
 
 
   Level: developer
 
 
 
.seealso: STSetBilinearForm(), STApplyB(), VecMDot(), STInnerProduct()
 
@*/
 
PetscErrorCode STMInnerProduct(ST st,PetscInt n,Vec x,const Vec y[],PetscScalar *p)
 
{
 
  PetscErrorCode ierr;
 
 
 
  PetscFunctionBegin;
 
  PetscValidHeaderSpecific(st,ST_COOKIE,1);
 
  PetscValidHeaderSpecific(x,VEC_COOKIE,3);
 
  PetscValidPointer(y,4);
 
  PetscValidHeaderSpecific(*y,VEC_COOKIE,4);
 
  PetscValidScalarPointer(p,5);
 
 
 
  ierr = PetscLogEventBegin(ST_InnerProduct,st,x,0,0);CHKERRQ(ierr);
 
  st->innerproducts += n;
 
  switch (st->bilinear_form) {
 
  case STINNER_HERMITIAN:
 
  case STINNER_SYMMETRIC:
 
    ierr = VecCopy(x,st->w);CHKERRQ(ierr);
 
    break;
 
  case STINNER_B_HERMITIAN:
 
  case STINNER_B_SYMMETRIC:
 
    ierr = STApplyB(st,x,st->w);CHKERRQ(ierr);
 
    break;
 
  }
 
  switch (st->bilinear_form) {
 
  case STINNER_HERMITIAN:
 
  case STINNER_B_HERMITIAN:
 
    ierr = VecMDot(st->w,n,y,p);CHKERRQ(ierr);
 
    break;
 
  case STINNER_SYMMETRIC:
 
  case STINNER_B_SYMMETRIC:
 
    ierr = VecMTDot(st->w,n,y,p);CHKERRQ(ierr);
 
    break;
 
  }
 
  ierr = PetscLogEventEnd(ST_InnerProduct,st,x,0,0);CHKERRQ(ierr);
 
  PetscFunctionReturn(0);
 
}
 
 
 
#undef __FUNCT__  
 
#define __FUNCT__ "STMInnerProductBegin"
 
/*@
 
   STMInnerProductBegin - Starts a split phase multiple inner product computation.
 
 
 
   Input Parameters:
 
+  st - the spectral transformation context
 
.  n  - number of vectors in y
 
.  x  - the first input vector
 
.  y  - array of vectors
 
-  p  - where the result will go
 
 
 
   Level: developer
 
 
 
   Notes:
 
   Each call to STMInnerProductBegin() should be paired with a call to STMInnerProductEnd().
 
 
 
.seealso: STMInnerProductEnd(), STMInnerProduct(), STNorm(), STNormBegin(),
 
          STNormEnd(), STInnerProduct()
 
 
 
@*/
 
PetscErrorCode STMInnerProductBegin(ST st,PetscInt n,Vec x,const Vec y[],PetscScalar *p)
 
{
 
  PetscErrorCode ierr;
 
 
 
  PetscFunctionBegin;
 
  PetscValidHeaderSpecific(st,ST_COOKIE,1);
 
  PetscValidHeaderSpecific(x,VEC_COOKIE,3);
 
  PetscValidPointer(y,4);
 
  PetscValidHeaderSpecific(*y,VEC_COOKIE,4);
 
  PetscValidScalarPointer(p,5);
 
 
 
  ierr = PetscLogEventBegin(ST_InnerProduct,st,x,0,0);CHKERRQ(ierr);
 
  st->innerproducts += n;
 
  switch (st->bilinear_form) {
 
  case STINNER_HERMITIAN:
 
  case STINNER_SYMMETRIC:
 
    ierr = VecCopy(x,st->w);CHKERRQ(ierr);
 
    break;
 
  case STINNER_B_HERMITIAN:
 
  case STINNER_B_SYMMETRIC:
 
    ierr = STApplyB(st,x,st->w);CHKERRQ(ierr);
 
    break;
 
  }
 
  switch (st->bilinear_form) {
 
  case STINNER_HERMITIAN:
 
  case STINNER_B_HERMITIAN:
 
    ierr = VecMDotBegin(st->w,n,y,p);CHKERRQ(ierr);
 
    break;
 
  case STINNER_SYMMETRIC:
 
  case STINNER_B_SYMMETRIC:
 
    ierr = VecMTDotBegin(st->w,n,y,p);CHKERRQ(ierr);
 
    break;
 
  }
 
  ierr = PetscLogEventEnd(ST_InnerProduct,st,x,0,0);CHKERRQ(ierr);
 
  PetscFunctionReturn(0);
 
}
 
 
 
#undef __FUNCT__  
 
#define __FUNCT__ "STMInnerProductEnd"
 
/*@
 
   STMInnerProductEnd - Ends a split phase multiple inner product computation.
 
 
 
   Input Parameters:
 
+  st - the spectral transformation context
 
.  n  - number of vectors in y
 
.  x  - the first input vector
 
-  y  - array of vectors
 
 
 
   Output Parameter:
 
.  p - result of the inner products
 
 
 
   Level: developer
 
 
 
   Notes:
 
   Each call to STMInnerProductBegin() should be paired with a call to STMInnerProductEnd().
 
 
 
.seealso: STMInnerProductBegin(), STMInnerProduct(), STNorm(), STNormBegin(),
 
          STNormEnd(), STInnerProduct()
 
 
 
@*/
 
PetscErrorCode STMInnerProductEnd(ST st,PetscInt n,Vec x,const Vec y[],PetscScalar *p)
 
{
 
  PetscErrorCode ierr;
 
 
 
  PetscFunctionBegin;
 
  PetscValidHeaderSpecific(st,ST_COOKIE,1);
 
  PetscValidHeaderSpecific(x,VEC_COOKIE,3);
 
  PetscValidPointer(y,4);
 
  PetscValidHeaderSpecific(*y,VEC_COOKIE,4);
 
  PetscValidScalarPointer(p,5);
 
 
 
  ierr = PetscLogEventBegin(ST_InnerProduct,st,x,0,0);CHKERRQ(ierr);
 
  switch (st->bilinear_form) {
 
  case STINNER_HERMITIAN:
 
  case STINNER_B_HERMITIAN:
 
    ierr = VecMDotEnd(st->w,n,y,p);CHKERRQ(ierr);
 
    break;
 
  case STINNER_SYMMETRIC:
 
  case STINNER_B_SYMMETRIC:
 
    ierr = VecMTDotEnd(st->w,n,y,p);CHKERRQ(ierr);
 
    break;
 
  }
 
  ierr = PetscLogEventEnd(ST_InnerProduct,st,x,0,0);CHKERRQ(ierr);
 
  PetscFunctionReturn(0);
  PetscFunctionReturn(0);
}
}
 
 
#undef __FUNCT__  
#undef __FUNCT__  
#define __FUNCT__ "STSetUp"
#define __FUNCT__ "STSetUp"