Actual source code: matreg.c

petsc-3.8.4 2018-03-24
Report Typos and Errors

  2: /*
  3:      Mechanism for register PETSc matrix types
  4: */
  5:  #include <petsc/private/matimpl.h>

  7: PetscBool MatRegisterAllCalled = PETSC_FALSE;

  9: /*
 10:    Contains the list of registered Mat routines
 11: */
 12: PetscFunctionList MatList = 0;

 14: /*@C
 15:    MatSetType - Builds matrix object for a particular matrix type

 17:    Collective on Mat

 19:    Input Parameters:
 20: +  mat      - the matrix object
 21: -  matype   - matrix type

 23:    Options Database Key:
 24: .  -mat_type  <method> - Sets the type; use -help for a list
 25:     of available methods (for instance, seqaij)

 27:    Notes:
 28:    See "${PETSC_DIR}/include/petscmat.h" for available methods

 30:   Level: intermediate

 32: .keywords: Mat, MatType, set, method

 34: .seealso: PCSetType(), VecSetType(), MatCreate(), MatType, Mat
 35: @*/
 36: PetscErrorCode  MatSetType(Mat mat, MatType matype)
 37: {
 38:   PetscErrorCode ierr,(*r)(Mat);
 39:   PetscBool      sametype,found;
 40:   MatBaseName    names = MatBaseNameList;


 45:   while (names) {
 46:     PetscStrcmp(matype,names->bname,&found);
 47:     if (found) {
 48:       PetscMPIInt size;
 49:       MPI_Comm_size(PetscObjectComm((PetscObject)mat),&size);
 50:       if (size == 1) matype = names->sname;
 51:       else matype = names->mname;
 52:       break;
 53:     }
 54:     names = names->next;
 55:   }

 57:   PetscObjectTypeCompare((PetscObject)mat,matype,&sametype);
 58:   if (sametype) return(0);

 60:    PetscFunctionListFind(MatList,matype,&r);
 61:   if (!r) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"Unknown Mat type given: %s",matype);

 63:   /* free the old data structure if it existed */
 64:   if (mat->ops->destroy) {
 65:     (*mat->ops->destroy)(mat);
 66:     mat->ops->destroy = NULL;

 68:     /* should these null spaces be removed? */
 69:     MatNullSpaceDestroy(&mat->nullsp);
 70:     MatNullSpaceDestroy(&mat->nearnullsp);
 71:     mat->preallocated = PETSC_FALSE;
 72:     mat->assembled = PETSC_FALSE;
 73:     mat->was_assembled = PETSC_FALSE;

 75:     /*
 76:      Increment, rather than reset these: the object is logically the same, so its logging and
 77:      state is inherited.  Furthermore, resetting makes it possible for the same state to be
 78:      obtained with a different structure, confusing the PC.
 79:     */
 80:     ++mat->nonzerostate;
 81:     PetscObjectStateIncrease((PetscObject)mat);
 82:   }
 83:   mat->preallocated  = PETSC_FALSE;
 84:   mat->assembled     = PETSC_FALSE;
 85:   mat->was_assembled = PETSC_FALSE;

 87:   /* increase the state so that any code holding the current state knows the matrix has been changed */
 88:   mat->nonzerostate++;
 89:   PetscObjectStateIncrease((PetscObject)mat);

 91:   /* create the new data structure */
 92:   (*r)(mat);
 93:   return(0);
 94: }

 96: /*@C
 97:    MatGetType - Gets the matrix type as a string from the matrix object.

 99:    Not Collective

101:    Input Parameter:
102: .  mat - the matrix

104:    Output Parameter:
105: .  name - name of matrix type

107:    Level: intermediate

109: .keywords: Mat, MatType, get, method, name

111: .seealso: MatSetType()
112: @*/
113: PetscErrorCode  MatGetType(Mat mat,MatType *type)
114: {
118:   *type = ((PetscObject)mat)->type_name;
119:   return(0);
120: }


123: /*@C
124:   MatRegister -  - Adds a new matrix type

126:    Not Collective

128:    Input Parameters:
129: +  name - name of a new user-defined matrix type
130: -  routine_create - routine to create method context

132:    Notes:
133:    MatRegister() may be called multiple times to add several user-defined solvers.

135:    Sample usage:
136: .vb
137:    MatRegister("my_mat",MyMatCreate);
138: .ve

140:    Then, your solver can be chosen with the procedural interface via
141: $     MatSetType(Mat,"my_mat")
142:    or at runtime via the option
143: $     -mat_type my_mat

145:    Level: advanced

147: .keywords: Mat, register

149: .seealso: MatRegisterAll()


152:   Level: advanced
153: @*/
154: PetscErrorCode  MatRegister(const char sname[],PetscErrorCode (*function)(Mat))
155: {

159:   PetscFunctionListAdd(&MatList,sname,function);
160:   return(0);
161: }

163: MatBaseName MatBaseNameList = 0;

165: /*@C
166:       MatRegisterBaseName - Registers a name that can be used for either a sequential or its corresponding parallel matrix type.

168:   Input Parameters:
169: +     bname - the basename, for example, MATAIJ
170: .     sname - the name of the sequential matrix type, for example, MATSEQAIJ
171: -     mname - the name of the parallel matrix type, for example, MATMPIAIJ


174:   Level: advanced
175: @*/
176: PetscErrorCode  MatRegisterBaseName(const char bname[],const char sname[],const char mname[])
177: {
179:   MatBaseName    names;

182:   PetscNew(&names);
183:   PetscStrallocpy(bname,&names->bname);
184:   PetscStrallocpy(sname,&names->sname);
185:   PetscStrallocpy(mname,&names->mname);
186:   if (!MatBaseNameList) {
187:     MatBaseNameList = names;
188:   } else {
189:     MatBaseName next = MatBaseNameList;
190:     while (next->next) next = next->next;
191:     next->next = names;
192:   }
193:   return(0);
194: }