Rocstar  1.0
Rocstar multiphysics simulation application
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
RFLU_ModPETScPoisson.F90
Go to the documentation of this file.
1 ! *********************************************************************
2 ! * Rocstar Simulation Suite *
3 ! * Copyright@2015, Illinois Rocstar LLC. All rights reserved. *
4 ! * *
5 ! * Illinois Rocstar LLC *
6 ! * Champaign, IL *
7 ! * www.illinoisrocstar.com *
8 ! * sales@illinoisrocstar.com *
9 ! * *
10 ! * License: See LICENSE file in top level of distribution package or *
11 ! * http://opensource.org/licenses/NCSA *
12 ! *********************************************************************
13 ! *********************************************************************
14 ! * THE SOFTWARE IS PROVIDED "AS IS", WITHOUT WARRANTY OF ANY KIND, *
15 ! * EXPRESS OR IMPLIED, INCLUDING BUT NOT LIMITED TO THE WARRANTIES *
16 ! * OF MERCHANTABILITY, FITNESS FOR A PARTICULAR PURPOSE AND *
17 ! * NONINFRINGEMENT. IN NO EVENT SHALL THE CONTRIBUTORS OR *
18 ! * COPYRIGHT HOLDERS BE LIABLE FOR ANY CLAIM, DAMAGES OR OTHER *
19 ! * LIABILITY, WHETHER IN AN ACTION OF CONTRACT, TORT OR OTHERWISE, *
20 ! * Arising FROM, OUT OF OR IN CONNECTION WITH THE SOFTWARE OR THE *
21 ! * USE OR OTHER DEALINGS WITH THE SOFTWARE. *
22 ! *********************************************************************
23 ! ******************************************************************************
24 !
25 ! Purpose: Suite of routines to solve Poisson equation using PETSc.
26 !
27 ! Description: None.
28 !
29 ! Notes: None.
30 !
31 ! ******************************************************************************
32 !
33 ! $Id: RFLU_ModPETScPoisson.F90,v 1.5 2008/12/06 08:44:23 mtcampbe Exp $
34 !
35 ! Copyright: (c) 2004 by the University of Illinois
36 !
37 ! ******************************************************************************
38 
40 
41  USE modglobal, ONLY: t_global
42  USE modparameters
43  USE moddatatypes
44  USE moderror
45  USE modbndpatch, ONLY: t_patch
46  USE moddatastruct, ONLY: t_region
47  USE modgrid, ONLY: t_grid
48  USE modmpi
49 
51 
52  IMPLICIT NONE
53 
54  PRIVATE
55  PUBLIC :: rflu_petsc_buildpoisson, &
59 
60 #include "include/finclude/petsc.h"
61 #include "include/finclude/petscvec.h"
62 #include "include/finclude/petscmat.h"
63 #include "include/finclude/petscksp.h"
64 #include "include/finclude/petscpc.h"
65 
66 ! ******************************************************************************
67 ! Declarations and definitions
68 ! ******************************************************************************
69 
70  CHARACTER(CHRLEN) :: &
71  RCSIdentString = '$RCSfile: RFLU_ModPETScPoisson.F90,v $ $Revision: 1.5 $'
72 
73 ! ******************************************************************************
74 ! Routines
75 ! ******************************************************************************
76 
77  CONTAINS
78 
79 
80 
81 
82 
83 
84 ! ******************************************************************************
85 !
86 ! Purpose: Build Poisson equation.
87 !
88 ! Description: None.
89 !
90 ! Input:
91 ! pRegion Pointer to region
92 !
93 ! Output: None.
94 !
95 ! Notes: None.
96 !
97 ! ******************************************************************************
98 
99  SUBROUTINE rflu_petsc_buildpoisson(pRegion)
100 
101  IMPLICIT NONE
102 
103 ! ******************************************************************************
104 ! Declarations and definitions
105 ! ******************************************************************************
106 
107 ! ==============================================================================
108 ! Parameters
109 ! ==============================================================================
110 
111  TYPE(t_region), POINTER :: pregion
112 
113 ! ==============================================================================
114 ! Locals
115 ! ==============================================================================
116 
117  ksp :: petsc_ksp
118  mat :: petsc_a
119  matnullspace :: petsc_nsp
120  pc :: petsc_pc
121  vec :: petsc_b,petsc_x
122  INTEGER :: petsc_col,petsc_row
123 
124  INTEGER :: c1,c2,errorflag,i,icg,ifg,ipatch,j,jbeg,jend
125  REAL(RFREAL) :: term
126  TYPE(t_global), POINTER :: global
127  TYPE(t_grid), POINTER :: pgrid
128  TYPE(t_patch), POINTER :: ppatch
129 
130 ! ******************************************************************************
131 ! Start
132 ! ******************************************************************************
133 
134  global => pregion%global
135 
136  CALL registerfunction(global,'RFLU_PETSC_BuildPoisson',&
137  'RFLU_ModPETScPoisson.F90')
138 
139  IF ( global%myProcid == masterproc .AND. &
140  global%verbLevel > verbose_low ) THEN
141  WRITE(stdout,'(A,1X,A)') solver_name,'Building Poisson...'
142  END IF ! global%verbLevel
143 
144  pgrid => pregion%grid
145 
146 ! ******************************************************************************
147 ! Initialize matrix
148 ! ******************************************************************************
149 
150  DO i = 1,SIZE(pgrid%poissonA,1)
151  pgrid%poissonA(i) = 0.0_rfreal
152  END DO ! i
153 
154  DO icg = 1,pgrid%nCells
155  pgrid%poissonLast(icg) = pgrid%poissonFirst(icg) + 1
156  END DO ! icg
157 
158 ! ******************************************************************************
159 ! Assemble matrix
160 ! ******************************************************************************
161 
162 ! ==============================================================================
163 ! Interior faces
164 ! ==============================================================================
165 
166  DO ifg = 1,pgrid%nFaces
167  c1 = pgrid%f2c(1,ifg)
168  c2 = pgrid%f2c(2,ifg)
169 
170  term = pgrid%fn(4,ifg)/(pgrid%cofgDist(1,ifg) + pgrid%cofgDist(2,ifg))
171 
172  IF ( c1 <= pgrid%nCells ) THEN
173  pgrid%poissonA(pgrid%poissonFirst(c1)) = &
174  pgrid%poissonA(pgrid%poissonFirst(c1)) + term
175  pgrid%poissonA(pgrid%poissonLast(c1)) = -term
176  pgrid%poissonLast(c1) = pgrid%poissonLast(c1) + 1
177  END IF ! c1
178 
179  IF ( c2 <= pgrid%nCells ) THEN
180  pgrid%poissonA(pgrid%poissonFirst(c2)) = &
181  pgrid%poissonA(pgrid%poissonFirst(c2)) + term
182  pgrid%poissonA(pgrid%poissonLast(c2)) = -term
183  pgrid%poissonLast(c2) = pgrid%poissonLast(c2) + 1
184  END IF ! c2
185  END DO ! ifg
186 
187 ! ==============================================================================
188 ! Boundary faces
189 ! ==============================================================================
190 
191 ! TO BE DONE
192 ! loop over Dirichlet boundary patches
193 ! but generally we do not have Dirichlet boundary conditions for pressure
194 ! END TO BE DONE
195 
196 ! ******************************************************************************
197 ! Set matrix. NOTE matrix lower boundary are always zero, even in FORTRAN.
198 ! ******************************************************************************
199 
200  petsc_a = pregion%poissonInfoPETSc(rflu_petsc_poisson_info_a)
201 
202  DO icg = 1,pgrid%nCells
203  jbeg = pgrid%poissonFirst(icg)
204  jend = pgrid%poissonFirst(icg) + pgrid%poissonNNZ(icg)-1
205 
206  DO j = jbeg,jend
207  petsc_row = icg - 1
208  petsc_col = pgrid%poissonCol(j) - 1
209 
210 
211  CALL matsetvalues(petsc_a,1,petsc_row,1,petsc_col, &
212  pgrid%poissonA(j),insert_values,errorflag)
213 
214  global%error = errorflag
215  IF ( global%error /= err_none ) THEN
216  CALL errorstop(global,err_petsc_output,__line__)
217  END IF ! global%error
218  END DO ! j
219  END DO ! icg
220 
221  CALL matassemblybegin(petsc_a,mat_final_assembly,errorflag)
222  global%error = errorflag
223  IF ( global%error /= err_none ) THEN
224  CALL errorstop(global,err_petsc_output,__line__)
225  END IF ! global%error
226 
227  CALL matassemblyend(petsc_a,mat_final_assembly,errorflag)
228  global%error = errorflag
229  IF ( global%error /= err_none ) THEN
230  CALL errorstop(global,err_petsc_output,__line__)
231  END IF ! global%error
232 
233 ! ******************************************************************************
234 ! End
235 ! ******************************************************************************
236 
237  IF ( global%myProcid == masterproc .AND. &
238  global%verbLevel > verbose_low ) THEN
239  WRITE(stdout,'(A,1X,A)') solver_name,'Building Poisson done.'
240  END IF ! global%verbLevel
241 
242  CALL deregisterfunction(global)
243 
244  END SUBROUTINE rflu_petsc_buildpoisson
245 
246 
247 
248 
249 
250 
251 
252 ! ******************************************************************************
253 !
254 ! Purpose:
255 !
256 ! Description: None.
257 !
258 ! Input:
259 ! pRegion Pointer to region
260 !
261 ! Output: None.
262 !
263 ! Notes: None.
264 !
265 ! ******************************************************************************
266 
267  SUBROUTINE rflu_petsc_createpoisson(pRegion)
268 
269  IMPLICIT NONE
270 
271 ! ******************************************************************************
272 ! Declarations and definitions
273 ! ******************************************************************************
274 
275 ! ==============================================================================
276 ! Parameters
277 ! ==============================================================================
278 
279  TYPE(t_region), POINTER :: pregion
280 
281 ! ==============================================================================
282 ! Locals
283 ! ==============================================================================
284 
285  mat :: petsc_a
286  vec :: petsc_b,petsc_x
287 
288  INTEGER :: c1,c2,errorflag,icg,ifg,ipatch
289  TYPE(t_global), POINTER :: global
290  TYPE(t_grid), POINTER :: pgrid
291  TYPE(t_patch), POINTER :: ppatch
292 
293 ! ******************************************************************************
294 ! Start
295 ! ******************************************************************************
296 
297  global => pregion%global
298 
299  CALL registerfunction(global,'RFLU_PETSC_CreatePoisson',&
300  'RFLU_ModPETScPoisson.F90')
301 
302  IF ( global%myProcid == masterproc .AND. &
303  global%verbLevel > verbose_low ) THEN
304  WRITE(stdout,'(A,1X,A)') solver_name,'Creating Poisson...'
305  END IF ! global%verbLevel
306 
307  pgrid => pregion%grid
308 
309 ! ******************************************************************************
310 ! Allocate memory for non-zero pattern and access arrays
311 ! ******************************************************************************
312 
313  ALLOCATE(pgrid%poissonNNZ(pgrid%nCells),stat=errorflag)
314  global%error = errorflag
315  IF ( global%error /= err_none ) THEN
316  CALL errorstop(global,err_allocate,__line__,'pGrid%poissonNNZ')
317  END IF ! global%error
318 
319  ALLOCATE(pgrid%poissonFirst(pgrid%nCells),stat=errorflag)
320  global%error = errorflag
321  IF ( global%error /= err_none ) THEN
322  CALL errorstop(global,err_allocate,__line__,'pGrid%poissonFirst')
323  END IF ! global%error
324 
325  ALLOCATE(pgrid%poissonLast(pgrid%nCells),stat=errorflag)
326  global%error = errorflag
327  IF ( global%error /= err_none ) THEN
328  CALL errorstop(global,err_allocate,__line__,'pGrid%poissonLast')
329  END IF ! global%error
330 
331 ! ******************************************************************************
332 ! Set non-zero pattern and access arrays
333 ! ******************************************************************************
334 
335  DO icg = 1,pgrid%nCells
336  pgrid%poissonNNZ(icg) = 1
337  END DO ! icg
338 
339  DO ifg = 1,pgrid%nFaces
340  c1 = pgrid%f2c(1,ifg)
341  c2 = pgrid%f2c(2,ifg)
342 
343  IF ( c1 <= pgrid%nCells ) THEN
344  pgrid%poissonNNZ(c1) = pgrid%poissonNNZ(c1) + 1
345  END IF ! c1
346 
347  IF ( c2 <= pgrid%nCells ) THEN
348  pgrid%poissonNNZ(c2) = pgrid%poissonNNZ(c2) + 1
349  END IF ! c2
350  END DO ! ifg
351 
352  pgrid%poissonFirst(1) = 1
353 
354  DO icg = 2,pgrid%nCells
355  pgrid%poissonFirst(icg) = pgrid%poissonFirst(icg-1) &
356  + pgrid%poissonNNZ(icg-1)
357  END DO ! icg
358 
359  ALLOCATE(pgrid%poissonA(sum(pgrid%poissonNNZ)),stat=errorflag)
360  global%error = errorflag
361  IF ( global%error /= err_none ) THEN
362  CALL errorstop(global,err_allocate,__line__,'pGrid%poissonA')
363  END IF ! global%error
364 
365  ALLOCATE(pgrid%poissonCol(sum(pgrid%poissonNNZ)),stat=errorflag)
366  global%error = errorflag
367  IF ( global%error /= err_none ) THEN
368  CALL errorstop(global,err_allocate,__line__,'pGrid%poissonCol')
369  END IF ! global%error
370 
371  DO icg = 1,pgrid%nCells
372  pgrid%poissonCol(pgrid%poissonFirst(icg)) = icg
373 
374  pgrid%poissonLast(icg) = pgrid%poissonFirst(icg) + 1
375  END DO ! icg
376 
377  DO ifg = 1,pgrid%nFaces
378  c1 = pgrid%f2c(1,ifg)
379  c2 = pgrid%f2c(2,ifg)
380 
381  IF ( c1 <= pgrid%nCells ) THEN
382  pgrid%poissonCol(pgrid%poissonLast(c1)) = c2
383  pgrid%poissonLast(c1) = pgrid%poissonLast(c1) + 1
384  END IF ! c1
385 
386  IF (c2 <= pgrid%nCells) THEN
387  pgrid%poissonCol(pgrid%poissonLast(c2)) = c1
388  pgrid%poissonLast(c2) = pgrid%poissonLast(c2) + 1
389  END IF ! c2
390  END DO ! ifg
391 
392 ! ******************************************************************************
393 ! Set matrix structure in PETSc
394 ! ******************************************************************************
395 
396  CALL matcreateseqaij(petsc_comm_self,pgrid%nCells,pgrid%nCells, &
397  petsc_null_integer,pgrid%poissonNNZ,petsc_a, &
398  errorflag)
399 
400  global%error = errorflag
401  IF ( global%error /= err_none ) THEN
402  CALL errorstop(global,err_petsc_output,__line__)
403  END IF ! global%error
404 
405  CALL veccreateseqwitharray(petsc_comm_self,pgrid%nCells,petsc_null_scalar, &
406  petsc_b,errorflag)
407 
408  global%error = errorflag
409  IF ( global%error /= err_none ) THEN
410  CALL errorstop(global,err_petsc_output,__line__)
411  END IF ! global%error
412 
413  CALL vecduplicate(petsc_b,petsc_x,errorflag)
414 
415  global%error = errorflag
416  IF ( global%error /= err_none ) THEN
417  CALL errorstop(global,err_petsc_output,__line__)
418  END IF ! global%error
419 
420  pregion%poissonInfoPETSc(rflu_petsc_poisson_info_a) = petsc_a
421  pregion%poissonInfoPETSc(rflu_petsc_poisson_info_b) = petsc_b
422  pregion%poissonInfoPETSc(rflu_petsc_poisson_info_x) = petsc_x
423 
424 ! ******************************************************************************
425 ! End
426 ! ******************************************************************************
427 
428  IF ( global%myProcid == masterproc .AND. &
429  global%verbLevel > verbose_low ) THEN
430  WRITE(stdout,'(A,1X,A)') solver_name,'Creating Poisson done.'
431  END IF ! global%verbLevel
432 
433  CALL deregisterfunction(global)
434 
435  END SUBROUTINE rflu_petsc_createpoisson
436 
437 
438 
439 
440 
441 
442 
443 
444 ! ******************************************************************************
445 !
446 ! Purpose: Set solver options for Poisson equation.
447 !
448 ! Description: None.
449 !
450 ! Input:
451 ! pRegion Pointer to region
452 !
453 ! Output: None.
454 !
455 ! Notes: None.
456 !
457 ! ******************************************************************************
458 
459  SUBROUTINE rflu_petsc_setsolverpoisson(pRegion)
460 
461  IMPLICIT NONE
462 
463 ! ******************************************************************************
464 ! Declarations and definitions
465 ! ******************************************************************************
466 
467 ! ==============================================================================
468 ! Parameters
469 ! ==============================================================================
470 
471  TYPE(t_region), POINTER :: pregion
472 
473 ! ==============================================================================
474 ! Locals
475 ! ==============================================================================
476 
477  ksp :: petsc_ksp
478  mat :: petsc_a
479  pc :: petsc_pc
480  matnullspace :: petsc_nsp
481 
482  INTEGER :: errorflag
483  TYPE(t_global), POINTER :: global
484 
485 ! ******************************************************************************
486 ! Start
487 ! ******************************************************************************
488 
489  global => pregion%global
490 
491  CALL registerfunction(global,'RFLU_PETSC_SetSolverPoisson',&
492  'RFLU_ModPETScPoisson.F90')
493 
494  IF ( global%myProcid == masterproc .AND. &
495  global%verbLevel > verbose_low ) THEN
496  WRITE(stdout,'(A,1X,A)') solver_name,'Setting solver Poisson...'
497  END IF ! global%verbLevel
498 
499 ! ******************************************************************************
500 ! Set solver and preconditioner options
501 ! ******************************************************************************
502 
503  petsc_a = pregion%poissonInfoPETSc(rflu_petsc_poisson_info_a)
504 
505  CALL kspcreate(petsc_comm_self,petsc_ksp,errorflag)
506  global%error = errorflag
507  IF ( global%error /= err_none ) THEN
508  CALL errorstop(global,err_petsc_output,__line__)
509  END IF ! global%error
510 
511  CALL kspsetoperators(petsc_ksp,petsc_a,petsc_a,same_nonzero_pattern, &
512  errorflag)
513  global%error = errorflag
514  IF ( global%error /= err_none ) THEN
515  CALL errorstop(global,err_petsc_output,__line__)
516  END IF ! global%error
517 
518  CALL kspsetfromoptions(petsc_ksp,errorflag)
519  global%error = errorflag
520  IF ( global%error /= err_none ) THEN
521  CALL errorstop(global,err_petsc_output,__line__)
522  END IF ! global%error
523 
524 ! TEMPORARY, needs to be sorted out
525 ! CALL MatNullSpaceCreate(PETSC_COMM_SELF,PETSC_TRUE,0,PETSC_NULL_VEC, &
526 ! PETSC_nsp,errorFlag)
527 ! END TEMPORARY
528 
529 ! Addition by hzhao
530  CALL matnullspacecreate(petsc_comm_self,petsc_true,0,0, &
531  petsc_nsp,errorflag)
532 ! End addition by hzhao
533 
534  global%error = errorflag
535  IF ( global%error /= err_none ) THEN
536  CALL errorstop(global,err_petsc_output,__line__)
537  END IF ! global%error
538 
539  CALL kspsetnullspace(petsc_ksp,petsc_nsp,errorflag)
540  global%error = errorflag
541  IF ( global%error /= err_none ) THEN
542  CALL errorstop(global,err_petsc_output,__line__)
543  END IF ! global%error
544 
545  CALL kspsettype(petsc_ksp,kspgmres,errorflag)
546  global%error = errorflag
547  IF ( global%error /= err_none ) THEN
548  CALL errorstop(global,err_petsc_output,__line__)
549  END IF ! global%error
550 
551  CALL kspgetpc(petsc_ksp,petsc_pc,errorflag)
552  global%error = errorflag
553  IF ( global%error /= err_none ) THEN
554  CALL errorstop(global,err_petsc_output,__line__)
555  END IF ! global%error
556 
557  CALL pcsettype(petsc_pc,pcilu,errorflag)
558  global%error = errorflag
559  IF ( global%error /= err_none ) THEN
560  CALL errorstop(global,err_petsc_output,__line__)
561  END IF ! global%error
562 
563  pregion%poissonInfoPETSc(rflu_petsc_poisson_info_ksp) = petsc_ksp
564  pregion%poissonInfoPETSc(rflu_petsc_poisson_info_pc ) = petsc_pc
565  pregion%poissonInfoPETSc(rflu_petsc_poisson_info_nsp) = petsc_nsp
566 
567 ! ******************************************************************************
568 ! End
569 ! ******************************************************************************
570 
571  IF ( global%myProcid == masterproc .AND. &
572  global%verbLevel > verbose_low ) THEN
573  WRITE(stdout,'(A,1X,A)') solver_name,'Setting solver Poisson done.'
574  END IF ! global%verbLevel
575 
576  CALL deregisterfunction(global)
577 
578  END SUBROUTINE rflu_petsc_setsolverpoisson
579 
580 
581 
582 
583 
584 
585 
586 
587 ! ******************************************************************************
588 !
589 ! Purpose: Solve Poisson equation.
590 !
591 ! Description: None.
592 !
593 ! Input:
594 ! pRegion Pointer to region
595 !
596 ! Output: None.
597 !
598 ! Notes: None.
599 !
600 ! ******************************************************************************
601 
603 
604  IMPLICIT NONE
605 
606 ! ******************************************************************************
607 ! Declarations and definitions
608 ! ******************************************************************************
609 
610 ! ==============================================================================
611 ! Parameters
612 ! ==============================================================================
613 
614  TYPE(t_region), POINTER :: pregion
615 
616 ! ==============================================================================
617 ! Locals
618 ! ==============================================================================
619 
620  ksp :: petsc_ksp
621  mat :: petsc_a
622  pc :: petsc_pc
623  vec :: petsc_b,petsc_x
624 
625  INTEGER :: c1,c2,errorflag,ifg,ifl,ipatch
626  REAL(RFREAL) :: flx
627  petscreal, DIMENSION(:), ALLOCATABLE :: b,x
628  TYPE(t_global), POINTER :: global
629  TYPE(t_grid), POINTER :: pgrid
630  TYPE(t_patch), POINTER :: ppatch
631 ! Addition by hzhao
632  INTEGER :: cvmixtpres
633 ! End addition by hzhao
634 ! ******************************************************************************
635 ! Start
636 ! ******************************************************************************
637 
638  global => pregion%global
639 
640  CALL registerfunction(global,'RFLU_PETSC_SolvePressurePoisson',&
641  'RFLU_ModPETScPoisson.F90')
642 
643  pgrid => pregion%grid
644 
645 ! ******************************************************************************
646 ! Allocate memory
647 ! ******************************************************************************
648 
649  ALLOCATE(x(pgrid%nCells),stat=errorflag)
650  global%error = errorflag
651  IF ( global%error /= err_none ) THEN
652  CALL errorstop(global,err_allocate,__line__,'x')
653  END IF ! global%error
654 
655  ALLOCATE(b(pgrid%nCells),stat=errorflag)
656  global%error = errorflag
657  IF ( global%error /= err_none ) THEN
658  CALL errorstop(global,err_allocate,__line__,'b')
659  END IF ! global%error
660 
661 ! ******************************************************************************
662 ! Set right-hand side
663 ! ******************************************************************************
664  DO ifg = 1,pgrid%nFaces
665  c1 = pgrid%f2c(1,ifg)
666  c2 = pgrid%f2c(2,ifg)
667 
668 ! Comment by hzhao
669 ! flx = pRegion%mixt%vfMixt(ifg)*pGrid%fn(XCOORD,ifg) &
670 ! + pRegion%mixt%vfMixt(ifg)*pGrid%fn(YCOORD,ifg) &
671 ! + pRegion%mixt%vfMixt(ifg)*pGrid%fn(ZCOORD,ifg)
672 ! End comment by hzhao
673 ! Addition by hzhao
674  flx = pregion%mixt%vfMixt(ifg)*pgrid%fn(4,ifg)
675 ! End addition by hzhao
676 
677  IF ( c1 <= pgrid%nCells ) THEN
678  b(c1) = b(c1) + flx
679  END IF ! c1
680 
681  IF ( c2 <= pgrid%nCells ) THEN
682  b(c2) = b(c2) - flx
683  END IF ! c2
684  END DO ! ifg
685 
686  DO ipatch = 1,pgrid%nPatches
687  ppatch => pregion%patches(ipatch)
688 
689  DO ifl = 1,ppatch%nBFaces
690  c1 = ppatch%bf2c(ifl)
691 
692 ! Comment by hzhao
693 ! flx = pPatch%vfMixt(ifl)*pPatch%fn(XCOORD,ifl) &
694 ! + pPatch%vfMixt(ifl)*pPatch%fn(YCOORD,ifl) &
695 ! + pPatch%vfMixt(ifl)*pPatch%fn(ZCOORD,ifl)
696 ! End comment by hzhao
697 ! Addition by hzhao
698  flx = ppatch%vfMixt(ifl)*ppatch%fn(4,ifl) &
699 ! End addition by hzhao
700 
701  b(c1) = b(c1) + flx
702  END DO ! ifl
703  END DO ! iPatch
704 
705 ! ******************************************************************************
706 ! Set left-hand side
707 ! ******************************************************************************
708 
709 ! TO DO
710 !
711 ! END TO DO
712 
713 ! ******************************************************************************
714 ! Solve
715 ! ******************************************************************************
716 
717  petsc_a = pregion%poissonInfoPETSc(rflu_petsc_poisson_info_a)
718  petsc_b = pregion%poissonInfoPETSc(rflu_petsc_poisson_info_b)
719  petsc_x = pregion%poissonInfoPETSc(rflu_petsc_poisson_info_x)
720 
721  petsc_ksp = pregion%poissonInfoPETSc(rflu_petsc_poisson_info_ksp)
722  petsc_pc = pregion%poissonInfoPETSc(rflu_petsc_poisson_info_pc )
723 
724  CALL vecplacearray(petsc_b,b,errorflag)
725  global%error = errorflag
726  IF ( global%error /= err_none ) THEN
727  CALL errorstop(global,err_petsc_output,__line__)
728  END IF ! global%error
729 
730  CALL vecplacearray(petsc_x,x,errorflag)
731  global%error = errorflag
732  IF ( global%error /= err_none ) THEN
733  CALL errorstop(global,err_petsc_output,__line__)
734  END IF ! global%error
735 
736 ! Addition by hzhao
737 ! Solvability condition
738  b = b - sum(b)/SIZE(b,1)
739 ! End addition by hzhao
740 
741  CALL kspsolve(petsc_ksp,petsc_b,petsc_x,errorflag)
742  global%error = errorflag
743  IF ( global%error /= err_none ) THEN
744  CALL errorstop(global,err_petsc_output,__line__)
745  END IF ! global%error
746 
747 ! ******************************************************************************
748 ! Store pressure
749 ! ******************************************************************************
750 ! Addition by hzhao
751  cvmixtpres = rflu_getcvloc(global, fluid_model_incomp, cv_mixt_pres)
752  pregion%mixt%cv(cvmixtpres,1:pgrid%nCells) = x
753 ! End addition by hzhao
754 
755 ! ******************************************************************************
756 ! Deallocate memory
757 ! ******************************************************************************
758 
759  DEALLOCATE(x,stat=errorflag)
760  global%error = errorflag
761  IF ( global%error /= err_none ) THEN
762  CALL errorstop(global,err_deallocate,__line__,'x')
763  END IF ! global%error
764 
765  DEALLOCATE(b,stat=errorflag)
766  global%error = errorflag
767  IF ( global%error /= err_none ) THEN
768  CALL errorstop(global,err_deallocate,__line__,'b')
769  END IF ! global%error
770 
771 ! ******************************************************************************
772 ! End
773 ! ******************************************************************************
774 
775  CALL deregisterfunction(global)
776 
777  END SUBROUTINE rflu_petsc_solvepressurepoisson
778 
779 
780 
781 
782 
783 
784 ! ******************************************************************************
785 ! End
786 ! ******************************************************************************
787 
788 END MODULE rflu_modpetscpoisson
789 
790 
791 ! ******************************************************************************
792 !
793 ! RCS Revision history:
794 !
795 ! $Log: RFLU_ModPETScPoisson.F90,v $
796 ! Revision 1.5 2008/12/06 08:44:23 mtcampbe
797 ! Updated license.
798 !
799 ! Revision 1.4 2008/11/19 22:17:34 mtcampbe
800 ! Added Illinois Open Source License/Copyright
801 !
802 ! Revision 1.3 2006/04/07 15:19:20 haselbac
803 ! Removed tabs
804 !
805 ! Revision 1.2 2005/01/13 21:39:38 haselbac
806 ! Incorporated Hongs additions for testing
807 !
808 ! Revision 1.1 2004/12/19 15:40:45 haselbac
809 ! Initial revision
810 !
811 ! ******************************************************************************
812 
813 
814 
815 
816 
817 
818 
819 
820 
821 
INTEGER function rflu_getcvloc(global, fluidModel, var)
Tfloat sum() const
Return the sum of all the pixel values in an image.
Definition: CImg.h:13022
subroutine, public rflu_petsc_solvepressurepoisson(pRegion)
subroutine, public rflu_petsc_buildpoisson(pRegion)
subroutine, public rflu_petsc_setsolverpoisson(pRegion)
subroutine registerfunction(global, funName, fileName)
Definition: ModError.F90:449
unsigned char b() const
Definition: Color.h:70
blockLoc i
Definition: read.cpp:79
void int int REAL * x
Definition: read.cpp:74
j indices j
Definition: Indexing.h:6
**********************************************************************Rocstar Simulation Suite Illinois Rocstar LLC All rights reserved ****Illinois Rocstar LLC IL **www illinoisrocstar com **sales illinoisrocstar com WITHOUT WARRANTY OF ANY **EXPRESS OR INCLUDING BUT NOT LIMITED TO THE WARRANTIES **OF FITNESS FOR A PARTICULAR PURPOSE AND **NONINFRINGEMENT IN NO EVENT SHALL THE CONTRIBUTORS OR **COPYRIGHT HOLDERS BE LIABLE FOR ANY DAMAGES OR OTHER WHETHER IN AN ACTION OF TORT OR **Arising OUT OF OR IN CONNECTION WITH THE SOFTWARE OR THE **USE OR OTHER DEALINGS WITH THE SOFTWARE **********************************************************************INTERFACE SUBROUTINE knode jend
subroutine errorstop(global, errorCode, errorLine, addMessage)
Definition: ModError.F90:483
**********************************************************************Rocstar Simulation Suite Illinois Rocstar LLC All rights reserved ****Illinois Rocstar LLC IL **www illinoisrocstar com **sales illinoisrocstar com WITHOUT WARRANTY OF ANY **EXPRESS OR INCLUDING BUT NOT LIMITED TO THE WARRANTIES **OF FITNESS FOR A PARTICULAR PURPOSE AND **NONINFRINGEMENT IN NO EVENT SHALL THE CONTRIBUTORS OR **COPYRIGHT HOLDERS BE LIABLE FOR ANY DAMAGES OR OTHER WHETHER IN AN ACTION OF TORT OR **Arising OUT OF OR IN CONNECTION WITH THE SOFTWARE OR THE **USE OR OTHER DEALINGS WITH THE SOFTWARE **********************************************************************INTERFACE SUBROUTINE knode jbeg
subroutine deregisterfunction(global)
Definition: ModError.F90:469
subroutine, public rflu_petsc_createpoisson(pRegion)